Type | Intent | Optional | Attributes | Name | ||
---|---|---|---|---|---|---|
integer(kind=idx_t), | intent(in) | :: | ne | The number of elements in the mesh. |
||
integer(kind=idx_t), | intent(in) | :: | nn | The number of nodes in the mesh. |
||
integer(kind=idx_t), | intent(in), | dimension(*) | :: | eptr | The pair of arrays storing the mesh as described in Section 5.6 of the manual. |
|
integer(kind=idx_t), | intent(in), | dimension(*) | :: | eind | The pair of arrays storing the mesh as described in Section 5.6 of the manual. |
|
integer(kind=idx_t), | intent(in) | :: | ncommon | The number of common nodes that two elements must have in order to put an edge between them in the dual graph. |
||
integer(kind=idx_t), | intent(in) | :: | numflag | Used to indicate which numbering scheme is used for |
||
type(c_ptr), | intent(out) | :: | xadj | These arrays store the adjacency structure of the generated dual graph. The format of the adjacency structure is described in Section 5.5 of the manual. |
||
type(c_ptr), | intent(out) | :: | adjncy | These arrays store the adjacency structure of the generated dual graph. The format of the adjacency structure is described in Section 5.5 of the manual. |
METIS_OK
- Indicates that the function returned normally.
METIS_ERROR_INPUT
- Indicates an input error.
METIS_ERROR_MEMORY
- Indicates that it could not allocate the required memory.
METIS_ERROR
- Indicates some other type of error.
This function is used to generate the dual graph of a mesh.
To use the returned arrays xadj
and adjncy
, these must be first converted from
a C pointer to a Fortran pointer using the subroutine c_f_pointer(cptr,fptr,shape)
that assigns the target of the C pointer cptr
to the Fortran pointer fptr
and
specifies its shape. The shape
is an integer rank-one array, storing the size ne+1
in case of the dual graph. The size of the new adjncy
array is stored in the
last element of xadj
when using C-style numbering. An example is shown below.
Memory for the returned arrays xadj
and adjncy
is allocated by METIS' API in C
using the standard malloc
function. It is the responsibility of the application to free
this memory by calling free
. Therefore, METIS provides the METIS_Free function that is a wrapper to
C's free
function.
The code below generates the nodal graph of the following mesh: image
use iso_c_binding, only : c_ptr, c_f_pointer use metis_interface, only: idx_t, METIS_MeshToDual, METIS_Free integer(kind=idx_t), parameter :: ne = 3, nn = 8, npel = 4 integer(kind=idx_t) :: ierr, eptr(ne+1), eind(ne*npel), numflag, ncommon type(c_ptr) :: xadj, adjncy integer(kind=idx_t), dimension(:), pointer :: fxadj => null(), fadjncy => null() numflag = 0 ! C-style numbering ncommon = 2 ! 2 common nodes per edge eptr = [0,4,8,12] eind = [0,1,2,3,1,4,5,2,4,6,7,5] ! note four nodes per element ierr = METIS_MeshToDual(ne,nn,eptr,eind,ncommon,numflag,xadj,adjncy) call c_f_pointer(xadj,fxadj,shape=[ne+1]) ! xadj is of the size ne+1 call c_f_pointer(adjncy,fadjncy,shape=[fxadj(ne+1)]) ! last value in xadj is the size of adjncy !... use values in fxadj and fadjncy ... call METIS_Free(xadj) call METIS_Free(adjncy) end