This module contains routines for translating global array indices into their local counterparts (that is, the indices into the local data structure holding a given process' chunk of a distributed array). The MCT domain decomposition descriptors GlobalMap and GlobalSegMap are both supported. Indices can be translated one-at-a-time using the GlobalToLocalIndex routine or many at once using the GlobalToLocalIndices routine.
This module also provides facilities for setting the local row and column indices for a SparseMatrix through the GlobalToLocalMatrix routines.
INTERFACE:
module m_GlobalToLocalUSES:
No external modules are used in the declaration section of this module. implicit none private ! exceptPUBLIC MEMBER FUNCTIONS:
public :: GlobalToLocalIndex ! Translate Global to Local index ! (i.e. recover local index for a ! point from its global index). public :: GlobalToLocalIndices ! Translate Global to Local indices ! (i.e. recover local starts/lengths ! of distributed data segments). public :: GlobalToLocalMatrix ! Re-indexing of row or column ! indices for a SparseMatrix interface GlobalToLocalIndices ; module procedure & GlobalSegMapToIndices_, & ! local arrays of starts/lengths GlobalSegMapToNavigator_, & ! return local indices as Navigator GlobalSegMapToIndexArr_ end interface interface GlobalToLocalIndex ; module procedure & GlobalSegMapToIndex_, & GlobalMapToIndex_ end interface interface GlobalToLocalMatrix ; module procedure & GlobalSegMapToLocalMatrix_ end interfaceSEE ALSO:
The MCT modules {\tt m\_GlobalMap} and {m\_GlobalSegMap} for more information regarding MCT's domain decomposition descriptors. The MCT module {\tt m\_SparseMatrix} for more information regarding the {\tt SparseMatrix} datatype.REVISION HISTORY:
2Feb01 - J.W. Larson <[email protected]> - initial prototype
GlobalSegMapToIndices_() takes a user-supplied GlobalSegMap data type GSMap, which desribes a decomposition on the input MPI communicator corresponding to the Fortran INTEGER handle comm to translate the global directory of segment locations into local indices for referencing the on-pe storage of the mapped distributed data.
N.B.: This routine returns two allocated arrays--start(:) and length(:)--which must be deallocated once the user no longer needs them. Failure to do this will create a memory leak.
INTERFACE:
subroutine GlobalSegMapToIndices_(GSMap, comm, start, length)USES:
use m_mpif90 use m_die, only : MP_perr_die, die, warn use m_GlobalSegMap, only : GlobalSegMap use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg use m_GlobalSegMap, only : GlobalSegMap_nlseg => nlseg implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! Output GlobalSegMap integer, intent(in) :: comm ! communicator handleOUTPUT PARAMETERS:
integer,dimension(:), pointer :: start ! local segment start indices integer,dimension(:), pointer :: length ! local segment sizesREVISION HISTORY:
2Feb01 - J.W. Larson <[email protected]> - initial version
This INTEGER query function takes a user-supplied GlobalSegMap data type GSMap, which desribes a decomposition on the input MPI communicator corresponding to the Fortran INTEGER handle comm, and the input global index value i_g, and returns a positive local index value if the datum i_g. If the datum i_g is not stored on the local process ID, a value of -1 is returned.
INTERFACE:
integer function GlobalSegMapToIndex_(GSMap, i_g, comm)USES:
use m_mpif90 use m_die, only : MP_perr_die, die, warn use m_GlobalSegMap, only : GlobalSegMap use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg use m_GlobalSegMap, only : GlobalSegMap_nlseg => nlseg implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! Output GlobalSegMap integer, intent(in) :: i_g ! global index integer, intent(in) :: comm ! communicator handleREVISION HISTORY:
2Feb01 - J.W. Larson <[email protected]> - initial version
Given a GlobalSegMap data type GSMap and MPI communicator corresponding to the Fortran INTEGER handle comm, convert an array of global index values i_global() to an array of local index values i_local(). If the datum i_global(j) is not stored on the local process ID, then i_local(j) will be set to -1/
INTERFACE:
subroutine GlobalSegMapToIndexArr_(GSMap, i_global, i_local, nindex, comm)USES:
use m_stdio use m_mpif90 use m_die, only : MP_perr_die, die, warn use m_GlobalSegMap, only : GlobalSegMap use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg use m_GlobalSegMap, only : GlobalSegMap_nlseg => nlseg implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! Output GlobalSegMap integer, intent(in) :: i_global(:) ! global index integer, intent(out) :: i_local(:) ! local index integer, intent(in) :: nindex ! size of i_global() integer, intent(in) :: comm ! communicator handleREVISION HISTORY:
12-apr-2006 R. Loy <[email protected]> - initial version
This INTEGER query function takes as its input a user-supplied GlobalMap data type GMap, which desribes a decomposition on the input MPI communicator corresponding to the Fortran INTEGER handle comm, and the input global index value i_g, and returns a positive local index value if the datum i_g. If the datum i_g is not stored on the local process ID, a value of -1 is returned.
INTERFACE:
integer function GlobalMapToIndex_(GMap, i_g, comm)USES:
use m_mpif90 use m_die, only : MP_perr_die, die, warn use m_GlobalMap, only : GlobalMap implicit noneINPUT PARAMETERS:
type(GlobalMap), intent(in) :: GMap ! Input GlobalMap integer, intent(in) :: i_g ! global index integer, intent(in) :: comm ! communicator handleREVISION HISTORY:
2Feb01 - J.W. Larson <[email protected]> - initial version
This routine takes as its input takes a user-supplied GlobalSegMap data type GSMap, which desribes a decomposition on the input MPI communicator corresponding to the Fortran INTEGER handle comm, and returns the local segment start index and length information for referencing the on-pe storage of the mapped distributed data. These data are returned in the form of the output Navigator argument Nav.
N.B.: This routine returns a Navigator variable Nav, which must be deallocated once the user no longer needs it. Failure to do this will create a memory leak.
INTERFACE:
subroutine GlobalSegMapToNavigator_(GSMap, comm, oNav)USES:
use m_mpif90 use m_die, only : MP_perr_die, die, warn use m_GlobalSegMap, only : GlobalSegMap use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg use m_GlobalSegMap, only : GlobalSegMap_nlseg => nlseg use m_Navigator, only : Navigator use m_Navigator, only : Navigator_init => init implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! Input GlobalSegMap integer, intent(in) :: comm ! communicator handleOUTPUT PARAMETERS:
type(Navigator), intent(out) :: oNav ! Output NavigatorREVISION HISTORY:
2Feb01 - J.W. Larson <[email protected]> - initial version
This routine takes as its input a user-supplied GlobalSegMap domain decomposition GSMap, which describes the decomposition of either the rows or columns of the input/output SparseMatrix argument sMat on the communicator associated with the INTEGER handle comm, and to translate the global row or column indices of sMat into their local counterparts. The choice of either row or column is governed by the value of the input CHARACTER argument RCFlag. One sets this variable to either 'ROW' or 'row' to specify row re-indexing (which are stored in sMat and retrieved by indexing the attribute lrow), and 'COLUMN' or 'column' to specify column re-indexing (which are stored in sMat and retrieved by indexing the SparseMatrix attribute lcol).
INTERFACE:
subroutine GlobalSegMapToLocalMatrix_(sMat, GSMap, RCFlag, comm)USES:
use m_stdio use m_die, only : die use m_SparseMatrix, only : SparseMatrix use m_SparseMatrix, only : SparseMatrix_indexIA => indexIA use m_SparseMatrix, only : SparseMatrix_lsize => lsize use m_GlobalSegMap, only : GlobalSegMap implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! Input GlobalSegMap character(len=*), intent(in) :: RCFlag ! 'row' or 'column' integer, intent(in) :: comm ! communicator handleINPUT/OUTPUT PARAMETERS:
type(SparseMatrix), intent(inout) :: sMatSEE ALSO:
The MCT module m_SparseMatrix for more information about the SparseMatrix type and its storage of global and local row-and column indices.REVISION HISTORY:
3May01 - J.W. Larson <[email protected]> - initial version, which is _extremely_ slow, but safe. This must be re-examined later.