Consider the problem of the 1-dimensional decomposition of an array across multiple processes. If each process owns only one contiguous segment, then the GlobalMap (see m_GlobalMap or details) is sufficient to describe the decomposition. If, however, each process owns multiple, non-adjacent segments of the array, a more sophisticated approach is needed. The GlobalSegMap data type allows one to describe a one-dimensional decomposition of an array with each process owning multiple, non-adjacent segments of the array.
In the current implementation of the GlobalSegMap, there is no
santity check to guarantee that
INTERFACE:
module m_GlobalSegMap implicit none private ! exceptPUBLIC MEMBER FUNCTIONS:
public :: GlobalSegMap ! The class data structure public :: init ! Create public :: clean ! Destroy public :: comp_id ! Return component ID number public :: gsize ! Return global vector size (excl. halos) public :: GlobalStorage ! Return total number of points in map, ! including halo points (if present). public :: ProcessStorage ! Return local storage on a given process. public :: OrderedPoints ! Return grid points of a given process in ! MCT-assumed order. public :: lsize ! Return local--that is, on-process--storage ! size (incl. halos) public :: ngseg ! Return global number of segments public :: nlseg ! Return local number of segments public :: max_nlseg ! Return max local number of segments public :: active_pes ! Return number of pes with at least 1 ! datum, and if requested, a list of them. public :: peLocs ! Given an input list of point indices, ! return its (unique) process ID. public :: haloed ! Is the input GlobalSegMap haloed? public :: rank ! Rank which process owns a datum public :: Sort ! compute index permutation to re-order ! GlobalSegMap%start, GlobalSegMap%length, ! and GlobalSegMap%pe_loc public :: Permute ! apply index permutation to re-order ! GlobalSegMap%start, GlobalSegMap%length, ! and GlobalSegMap%pe_loc public :: SortPermute ! compute index permutation and apply it to ! re-order the GlobalSegMap components ! GlobalSegMap%start, GlobalSegMap%length, ! and GlobalSegMap%pe_loc public :: increasing ! Are the indices for each pe strictly ! increasing? public :: copy ! Copy the gsmap public :: print ! Print the contents of the GSMapPUBLIC TYPES:
type GlobalSegMap #ifdef SEQUENCE sequence #endif integer :: comp_id ! Component ID number integer :: ngseg ! No. of Global segments integer :: gsize ! No. of Global elements integer,dimension(:),pointer :: start ! global seg. start index integer,dimension(:),pointer :: length ! segment lengths integer,dimension(:),pointer :: pe_loc ! PE locations end type GlobalSegMap interface init ; module procedure & initd_, & ! initialize from all PEs initr_, & ! initialize from the root initp_, & ! initialize in parallel from replicated arrays initp1_, & ! initialize in parallel from 1 replicated array initp0_, & ! null constructor using replicated data init_index_ ! initialize from local index arrays end interface interface clean ; module procedure clean_ ; end interface interface comp_id ; module procedure comp_id_ ; end interface interface gsize ; module procedure gsize_ ; end interface interface GlobalStorage ; module procedure & GlobalStorage_ end interface interface ProcessStorage ; module procedure & ProcessStorage_ end interface interface OrderedPoints ; module procedure & OrderedPoints_ end interface interface lsize ; module procedure lsize_ ; end interface interface ngseg ; module procedure ngseg_ ; end interface interface nlseg ; module procedure nlseg_ ; end interface interface max_nlseg ; module procedure max_nlseg_ ; end interface interface active_pes ; module procedure active_pes_ ; end interface interface peLocs ; module procedure peLocs_ ; end interface interface haloed ; module procedure haloed_ ; end interface interface rank ; module procedure & rank1_ , & ! single rank case rankm_ ! degenerate (multiple) ranks for halo case end interface interface Sort ; module procedure Sort_ ; end interface interface Permute ; module procedure & PermuteInPlace_ end interface interface SortPermute ; module procedure & SortPermuteInPlace_ end interface interface increasing ; module procedure increasing_ ; end interface interface copy ; module procedure copy_ ; end interface interface print ; module procedure & print_ ,& printFromRootnp_ end interfaceREVISION HISTORY:
28Sep00 - J.W. Larson <[email protected]> - initial prototype 26Jan01 - J.W. Larson <[email protected]> - replaced the component GlobalSegMap%comm with GlobalSegMap%comp_id. 06Feb01 - J.W. Larson <[email protected]> - removed the GlobalSegMap%lsize component. Also, added the GlobalStorage query function. 24Feb01 - J.W. Larson <[email protected]> - Added the replicated initialization routines initp_() and initp1(). 25Feb01 - J.W. Larson <[email protected]> - Added the routine ProcessStorage_(). 18Apr01 - J.W. Larson <[email protected]> - Added the routine peLocs(). 26Apr01 - R. Jacob <[email protected]> - Added the routine OrderedPoints_(). 03Aug01 - E. Ong <[email protected]> - In initd_, call initr_ with actual shaped arguments on non-root processes to satisfy F90 standard. See comments in initd. 18Oct01 - J.W. Larson <[email protected]> - Added the routine bcast(), and also cleaned up prologues.
This routine takes the scattered input INTEGER arrays start, length, and pe_loc, gathers these data to the root process, and from them creates a global set of segment information for the output GlobalSegMap argument GSMap. The input INTEGER arguments comp_id, gsize provide the GlobalSegMap component ID number and global grid size, respectively. The input argument my_comm is the F90 INTEGER handle for the MPI communicator. If the input arrays are overdimensioned, optional argument numel can be used to specify how many elements should be used.
INTERFACE:
subroutine initd_(GSMap, start, length, root, my_comm, & comp_id, pe_loc, gsize, numel)USES:
use m_mpif90 use m_die use m_stdio use m_FcComms, only : fc_gather_int, fc_gatherv_int implicit noneINPUT PARAMETERS:
integer,dimension(:),intent(in) :: start ! segment local start ! indices integer,dimension(:),intent(in) :: length ! segment local lengths integer,intent(in) :: root ! root on my_com integer,intent(in) :: my_comm ! local communicatior integer,intent(in) :: comp_id ! component model ID integer,dimension(:), pointer, optional :: pe_loc ! process location integer,intent(in), optional :: gsize ! global vector size ! (optional). It can ! be computed by this ! routine if no haloing ! is assumed. integer,intent(in), optional :: numel ! specify number of elements ! to use in start, lengthOUTPUT PARAMETERS:
type(GlobalSegMap),intent(out) :: GSMap ! Output GlobalSegMapREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype 14Nov00 - J.W. Larson <[email protected]> - final working version 09Jan01 - J.W. Larson <[email protected]> - repaired: a subtle bug concerning the usage of the argument pe_loc (result was the new pointer variable my_pe_loc); a mistake in the tag arguments to MPI_IRECV; a bug in the declaration of the array status used by MPI_WAITALL. 26Jan01 - J.W. Larson <[email protected]> - replaced optional argument gsm_comm with required argument comp_id. 23Sep02 - Add optional argument numel to allow start, length arrays to be overdimensioned. 31Jan09 - P.H. Worley <[email protected]> - replaced irecv/send/waitall logic with calls to flow controlled gather routines
This routine takes the input INTEGER arrays start, length, and pe_loc (all valid only on the root process), and from them creates a global set of segment information for the output GlobalSegMap argument GSMap. The input INTEGER arguments ngseg, comp_id, gsize (again, valid only on the root process) provide the GlobalSegMap global segment count, component ID number, and global grid size, respectively. The input argument my_comm is the F90 INTEGER handle for the MPI communicator.
INTERFACE:
subroutine initr_(GSMap, ngseg, start, length, pe_loc, root, & my_comm, comp_id, gsize)USES:
use m_mpif90 use m_die use m_stdio implicit noneINPUT PARAMETERS:
integer, intent(in) :: ngseg ! no. of global segments integer,dimension(:),intent(in) :: start ! segment local start index integer,dimension(:),intent(in) :: length ! the distributed sizes integer,dimension(:),intent(in) :: pe_loc ! process location integer,intent(in) :: root ! root on my_com integer,intent(in) :: my_comm ! local communicatior integer,intent(in) :: comp_id ! component id number integer,intent(in), optional :: gsize ! global vector size ! (optional). It can ! be computed by this ! routine if no haloing ! is assumed.OUTPUT PARAMETERS:
type(GlobalSegMap),intent(out) :: GSMap ! Output GlobalSegMapREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype 09Nov00 - J.W. Larson <[email protected]> - final working version 10Jan01 - J.W. Larson <[email protected]> - minor bug fix 12Jan01 - J.W. Larson <[email protected]> - minor bug fix regarding disparities in ngseg on the root and other processes 26Jan01 - J.W. Larson <[email protected]> - replaced optional argument gsm_comm with required argument comp_id.
The routine initp_() takes the input replicated arguments comp_id, ngseg, gsize, start(:), length(:), and pe_loc(:), and uses them to initialize an output GlobalSegMap GSMap. This routine operates on the assumption that these data are replicated across the communicator on which the GlobalSegMap is being created.
INTERFACE:
subroutine initp_(GSMap, comp_id, ngseg, gsize, start, length, pe_loc)USES:
use m_mpif90 use m_die, only : die use m_stdio implicit noneINPUT PARAMETERS:
integer,intent(in) :: comp_id ! component model ID integer,intent(in) :: ngseg ! global number of segments integer,intent(in) :: gsize ! global vector size integer,dimension(:),intent(in) :: start ! segment local start index integer,dimension(:),intent(in) :: length ! the distributed sizes integer,dimension(:),intent(in) :: pe_loc ! process locationOUTPUT PARAMETERS:
type(GlobalSegMap),intent(out) :: GSMap ! Output GlobalSegMapREVISION HISTORY:
24Feb01 - J.W. Larson <[email protected]> - Initial version.
The routine initp1_() takes the input replicated arguments
comp_id, ngseg, gsize, and all_arrays(:),
and uses them to initialize an output GlobalSegMap GSMap.
This routine operates on the assumption that these data are replicated
across the communicator on which the GlobalSegMap is being created.
The input array all_arrays(:) should be of length 2 * ngseg,
and is packed so that
INTERFACE:
subroutine initp1_(GSMap, comp_id, ngseg, gsize, all_arrays)USES:
use m_mpif90 use m_die, only : die use m_stdio implicit noneINPUT PARAMETERS:
integer,intent(in) :: comp_id ! component model ID integer,intent(in) :: ngseg ! global no. of segments integer,intent(in) :: gsize ! global vector size integer,dimension(:),intent(in) :: all_arrays ! packed array of length ! 3*ngseg containing (in ! this order): start(:), ! length(:), and pe_loc(:)OUTPUT PARAMETERS:
type(GlobalSegMap),intent(out) :: GSMap ! Output GlobalSegMapREVISION HISTORY:
24Feb01 - J.W. Larson <[email protected]> - Initial version.
The routine initp0_() takes the input replicated arguments comp_id, ngseg, gsize, and uses them perform null construction of the output GlobalSegMap GSMap. This is a null constructor in the sense that we are not filling in the segment information arrays. This routine operates on the assumption that these data are replicated across the communicator on which the GlobalSegMap is being created.
INTERFACE:
subroutine initp0_(GSMap, comp_id, ngseg, gsize)USES:
use m_die, only : die use m_stdio implicit noneINPUT PARAMETERS:
integer,intent(in) :: comp_id ! component model ID integer,intent(in) :: ngseg ! global number of segments integer,intent(in) :: gsize ! global vector sizeOUTPUT PARAMETERS:
type(GlobalSegMap),intent(out) :: GSMap ! Output GlobalSegMapREVISION HISTORY:
13Aug03 - J.W. Larson <[email protected]> - Initial version.
The routine init_index_() takes a local array of indices lindx and uses them to create a GlobalSegMap. lindx is parsed to determine the lengths of the runs, and then a call is made to initd_. The optional argument lsize can be used if only the first lsize number of elements of lindx are valid. The optional argument gsize is used to specify the global number of unique points if this can not be determined from the collective lindx.
INTERFACE:
subroutine init_index_(GSMap, lindx, my_comm, comp_id, lsize, gsize)USES:
use m_GlobalSegMap,only: GlobalSegMap use m_GlobalSegMap,only: MCT_GSMap_init => init use shr_sys_mod use m_die implicit noneINPUT PARAMETERS:
integer , dimension(:),intent(in) :: lindx ! index buffer integer , intent(in) :: my_comm ! mpi communicator group (mine) integer , intent(in) :: comp_id ! component id (mine) integer , intent(in),optional :: lsize ! size of index buffer integer , intent(in),optional :: gsize ! global vector sizeOUTPUT PARAMETERS:
type(GlobalSegMap),intent(out) :: GSMap ! Output GlobalSegMapREVISION HISTORY:
30Jul02 - T. Craig - initial version in cpl6. 17Nov05 - R. Loy <[email protected]> - install into MCT 18Nov05 - R. Loy <[email protected]> - make lsize optional 25Jul06 - R. Loy <[email protected]> - error check on lindex/alloc/dealloc
This routine deallocates the array components of the GlobalSegMap argument GSMap: GSMap%start, GSMap%length, and GSMap%pe_loc. It also zeroes out the values of the integer components GSMap%ngseg, GSMap%comp_id, and GSMap%gsize.
INTERFACE:
subroutine clean_(GSMap,stat)USES:
use m_die implicit noneINPUT/OUTPUT PARAMETERS:
type(GlobalSegMap), intent(inout) :: GSMap integer, optional, intent(out) :: statREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype 01Mar02 - E.T. Ong <[email protected]> - added stat argument. Removed dies to prevent crashing.
The function ngseg_() returns the global number of vector segments in the GlobalSegMap argument GSMap. This is merely the value of GSMap%ngseg.
INTERFACE:
integer function ngseg_(GSMap) implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: GSMapREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype
The function nlseg_() returns the number of vector segments in the GlobalSegMap argument GSMap that reside on the process specified by the input argument pID. This is the number of entries GSMap%pe_loc whose value equals pID.
INTERFACE:
integer function nlseg_(GSMap, pID) implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: GSMap integer, intent(in) :: pIDREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype 14Jun01 - J.W. Larson <[email protected]> - Bug fix in lower limit of loop over elements of GSMap%pe_loc(:). The original code had this lower limit set to 0, which was out-of-bounds (but uncaught). The correct lower index is 1. This bug was discovered by Everest Ong.
The function max_nlseg_() returns the maximum number over all processors of the vector segments in the GlobalSegMap argument gsap E.g. max_p(nlseg(gsmap,p)) but computed more efficiently
INTERFACE:
integer function max_nlseg_(gsmap)USES:
use m_MCTWorld, only :ThisMCTWorld use m_mpif90 use m_die use m_stdio ! rml implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: gsmapREVISION HISTORY:
17Jan07 - R. Loy <[email protected]> - initial prototype
The function comp_id_() returns component ID number stored in GSMap%comp_id.
INTERFACE:
integer function comp_id_(GSMap)USES:
use m_die,only: die use m_stdio, only :stderr implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: GSMapREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype 26Jan01 - J.W. Larson <[email protected]> - renamed comp_id_ to fit within MCT_World component ID context. 01May01 - R.L. Jacob <[email protected]> - make sure GSMap is defined.
The function gsize_() takes the input GlobalSegMap arguement GSMap and returns the global vector length stored in GlobalSegMap%gsize.
INTERFACE:
integer function gsize_(GSMap) implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: GSMapREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype
The function GlobalStorage_() takes the input GlobalSegMap arguement GSMap and returns the global storage space required (i.e., the vector length) to hold all the data specified by GSMap.
N.B.: If GSMap contains halo or masked points, the value by GlobalStorage_() may differ from GSMap%gsize.
INTERFACE:
integer function GlobalStorage_(GSMap) implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: GSMapREVISION HISTORY:
06Feb01 - J.W. Larson <[email protected]> - initial version
The function ProcessStorage_() takes the input GlobalSegMap arguement GSMap and returns the storage space required by process PEno (i.e., the vector length) to hold all the data specified by GSMap.
INTERFACE:
integer function ProcessStorage_(GSMap, PEno) implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: GSMap integer, intent(in) :: PEnoREVISION HISTORY:
06Feb01 - J.W. Larson <[email protected]> - initial version
returned in the assumed MCT order.
The function OrderedPoints_() takes the input GlobalSegMap arguement GSMap and returns a vector of the points owned by PEno. Points is allocated here. The calling process is responsible for deallocating the space.
INTERFACE:
subroutine OrderedPoints_(GSMap, PEno, Points)USES:
use m_die,only: die implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! input GlobalSegMap integer, intent(in) :: PEno ! input process number integer,dimension(:),pointer :: Points ! the vector of pointsREVISION HISTORY:
25Apr01 - R. Jacob <[email protected]> - initial prototype
This function returns the number of points owned by the local process, as defined by the input GlobalSegMap argument GSMap. The local process ID is determined through use of the input INTEGER argument comm, which is the Fortran handle for the MPI communicator.
INTERFACE:
integer function lsize_(GSMap, comm)USES:
use m_mpif90 use m_die , only : MP_perr_die implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap integer, intent(in) :: commREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype 06Feb01 - J.W. Larson <[email protected]> - Computed directly from the GlobalSegMap, rather than returning a hard- wired local attribute. This required the addition of the communicator argument.
index.
This routine assumes that there is one process that owns the datum with a given global index. It should not be used when the input GlobalSegMap argument GSMap has been built to incorporate halo points.
INTERFACE:
subroutine rank1_(GSMap, i_g, rank) implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! input GlobalSegMap integer, intent(in) :: i_g ! a global indexOUTPUT PARAMETERS:
integer, intent(out) :: rank ! the pe on which this ! element residesREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype
index.
This routine assumes that there may be more than one process that owns the datum with a given global index. This routine should be used when the input GlobalSegMap argument GSMap has been built to incorporate ! halo points. Nota Bene: The output array rank is allocated in this routine and must be deallocated by the routine calling rankm_(). Failure to do so could result in a memory leak.
INTERFACE:
subroutine rankm_(GSMap, i_g, num_loc, rank) implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! input GlobalSegMap integer, intent(in) :: i_g ! a global indexOUTPUT PARAMETERS:
integer, intent(out) :: num_loc ! the number of processes ! which own element i_g integer, dimension(:), pointer :: rank ! the process(es) on which ! element i_g residesREVISION HISTORY:
29Sep00 - J.W. Larson <[email protected]> - initial prototype
index.
This routine scans the pe location list of the input GlobalSegMap GSMap%pe_loc(:), and counts the number of pe locations that own at least one datum. This value is returned in the INTEGER argument n_active. If the optional INTEGER array argument list is included in the call, a sorted list (in ascending order) of the active processes will be returned.
N.B.: If active_pes_() is invoked with the optional argument pe_list included, this routine will allocate and return this array. The user must deallocate this array once it is no longer needed. Failure to do so will result in a memory leak.
INTERFACE:
subroutine active_pes_(GSMap, n_active, pe_list)USES:
use m_die , only : die use m_SortingTools , only : IndexSet use m_SortingTools , only : IndexSort use m_SortingTools , only : Permute implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMapOUTPUT PARAMETERS:
integer, intent(out) :: n_active integer, dimension(:), pointer, optional :: pe_listREVISION HISTORY:
03Feb01 - J.W. Larson <[email protected]> - initial version.
index.
This routine takes an input INTEGER array of point indices points(:), compares them with an input GlobalSegMap pointGSMap, and returns the unique process ID location for each point. Note the emphasize on unique. The assumption here (which is tested) is that pointGSMap is not haloed. The process ID locations for the points is returned in the array pe_locs(:).
N.B.: The test of pointGSMap for halo points, and the subsequent search for the process ID for each point is very slow. This first version of the routine is serial. A parallel version of this routine will need to be developed.
INTERFACE:
subroutine peLocs_(pointGSMap, npoints, points, pe_locs)USES:
use m_die , only : die implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: pointGSMap integer, intent(in) :: npoints integer, dimension(:), intent(in) :: pointsOUTPUT PARAMETERS:
integer, dimension(:), intent(out) :: pe_locsREVISION HISTORY:
18Apr01 - J.W. Larson <[email protected]> - initial version.
index.
This LOGICAL function tests the input GlobalSegMap GSMap for the presence of halo points. Halo points are points that appear in more than one segment of a GlobalSegMap. If any halo point is found, the function haloed_() returns immediately with value .TRUE. If, after an exhaustive search of the map has been completed, no halo points are found, the function haloed_() returns with value .FALSE.
The search algorithm is:
N.B.: Beware that the search for halo points is potentially expensive.
INTERFACE:
logical function haloed_(GSMap)USES:
use m_die , only : die use m_SortingTools , only : IndexSet use m_SortingTools , only : IndexSort use m_SortingTools , only : Permute implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMapREVISION HISTORY:
08Feb01 - J.W. Larson <[email protected]> - initial version. 26Apr01 - J.W. Larson <[email protected]> - Bug fix.
Sort_() uses the supplied keys key1 and key2 to generate a permutation perm that will put the entries of the components GlobalSegMap%start, GlobalSegMap%length and GlobalSegMap%pe_loc in ascending lexicographic order.
N.B.: Sort_() returns an allocated array perm(:). It the user must deallocate this array once it is no longer needed. Failure to do so could create a memory leak.
INTERFACE:
subroutine Sort_(GSMap, key1, key2, perm)USES:
use m_die , only : die use m_SortingTools , only : IndexSet use m_SortingTools , only : IndexSort implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: GSMap ! input GlobalSegMap integer, dimension(:), intent(in) :: key1 ! first sort key integer, dimension(:), intent(in), optional :: key2 ! second sort keyOUTPUT PARAMETERS:
integer, dimension(:), pointer :: perm ! output index permutationREVISION HISTORY:
02Feb01 - J.W. Larson <[email protected]> - initial version
PermuteInPlace_() uses a supplied index permutation perm to re-order GlobalSegMap%start, GlobalSegMap%length and GlobalSegMap%pe_loc.
INTERFACE:
subroutine PermuteInPlace_(GSMap, perm)USES:
use m_die , only : die use m_SortingTools , only : Permute implicit noneINPUT PARAMETERS:
integer, dimension(:), intent(in) :: permINPUT/OUTPUT PARAMETERS:
type(GlobalSegMap), intent(inout) :: GSMapREVISION HISTORY:
02Feb01 - J.W. Larson <[email protected]> - initial version.
SortPermuteInPlace_() uses a the supplied key(s) to generate and apply an index permutation that will place the GlobalSegMap components GlobalSegMap%start, GlobalSegMap%length and GlobalSegMap%pe_loc in lexicographic order.
INTERFACE:
subroutine SortPermuteInPlace_(GSMap, key1, key2)USES:
use m_die , only : die implicit noneINPUT PARAMETERS:
integer, dimension(:), intent(in) :: key1 integer, dimension(:), intent(in), optional :: key2INPUT/OUTPUT PARAMETERS:
type(GlobalSegMap), intent(inout) :: GSMapREVISION HISTORY:
02Feb01 - J.W. Larson <[email protected]> - initial version.
The function increasing_() returns .TRUE. if each proc's indices in the GlobalSegMap argument GSMap have strictly increasing indices. I.e. the proc's segments have indices in ascending order and are non-overlapping.
INTERFACE:
logical function increasing_(gsmap)USES:
use m_MCTWorld, only: ThisMCTWorld use m_die implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: gsmapREVISION HISTORY:
06Jun07 - R. Loy <[email protected]> - initial version
Make a copy of a gsmap. Note this is a deep copy of all arrays.
INTERFACE:
subroutine copy_(src,dest)USES:
use m_MCTWorld, only: ThisMCTWorld use m_die implicit noneINPUT PARAMETERS:
type(GlobalSegMap),intent(in) :: srcOUTPUT PARAMETERS:
type(GlobalSegMap),intent(out) :: destREVISION HISTORY:
27Jul07 - R. Loy <[email protected]> - initial version
Print out contents of GSMAP on unit number 'lun'
INTERFACE:
subroutine print_(gsmap,lun)USES:
use m_die implicit noneINPUT/OUTPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: gsmap integer, intent(in) :: lunREVISION HISTORY:
06Jul12 - R. Jacob <[email protected]> - initial version
Print out contents of GSMAP on unit number 'lun'
INTERFACE:
subroutine printFromRootnp_(gsmap,mycomm,lun)USES:
use m_MCTWorld, only : printnp use m_die use m_mpif90 implicit noneINPUT/OUTPUT PARAMETERS:
type(GlobalSegMap), intent(in) :: gsmap integer, intent(in) :: mycomm integer, intent(in) :: lunREVISION HISTORY:
06Jul12 - R. Jacob <[email protected]> - initial version
This module provides communications support for the GlobalSegMap datatype. Both blocking and non-blocking point-to-point communications are provided for send (analogues to MPI_SEND()/MPI_ISEND()) A receive and broadcast method is also supplied.
INTERFACE:
module m_GlobalSegMapComms implicit none private ! exceptPUBLIC MEMBER FUNCTIONS:
public :: send public :: recv public :: isend public :: bcast interface bcast ; module procedure bcast_ ; end interface interface send ; module procedure send_ ; end interface interface recv ; module procedure recv_ ; end interface interface isend ; module procedure isend_ ; end interfaceREVISION HISTORY:
11Aug03 - J.W. Larson <[email protected]> - initial version
This routine performs a blocking send of a GlobalSegMap (the input argument outgoingGSMap) to the root processor on component comp_id. The input INTEGER argument TagBase is used to generate tags for the messages associated with this operation; there are six messages involved, so the user should avoid using tag values TagBase and TagBase + 5. All six messages are blocking. The success (failure) of this operation is reported in the zero (non-zero) value of the optional INTEGER output variable status.
INTERFACE:
subroutine send_(outgoingGSMap, comp_id, TagBase, status)USES:
use m_mpif90 use m_die, only : MP_perr_die,die use m_stdio use m_GlobalSegMap, only : GlobalSegMap use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_ID use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize use m_MCTWorld, only : ComponentToWorldRank use m_MCTWorld, only : ThisMCTWorld implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(IN) :: outgoingGSMap integer, intent(IN) :: comp_id integer, intent(IN) :: TagBaseOUTPUT PARAMETERS:
integer, optional, intent(OUT) :: statusREVISION HISTORY:
13Aug03 - J.W. Larson <[email protected]> - API and initial version. 26Aug03 - R. Jacob <[email protected]> - use same method as isend_ 05Mar04 - R. Jacob <[email protected]> - match new isend_ method.
This routine performs a non-blocking send of a GlobalSegMap (the input argument outgoingGSMap) to the root processor on component comp_id The input INTEGER argument TagBase is used to generate tags for the messages associated with this operation; there are six messages involved, so the user should avoid using tag values TagBase and TagBase + 5. All six messages are non- blocking, and the request handles for them are returned in the output INTEGER array reqHandle, which can be checked for completion using any of MPI's wait functions. The success (failure) of this operation is reported in the zero (non-zero) value of the optional INTEGER output variable status.
N.B.: Data is sent directly out of outgoingGSMap so it must not be deleted until the send has completed.
N.B.: The array reqHandle represents allocated memory that must be deallocated when it is no longer needed. Failure to do so will create a memory leak.
INTERFACE:
subroutine isend_(outgoingGSMap, comp_id, TagBase, reqHandle, status)USES:
use m_mpif90 use m_die, only : MP_perr_die,die use m_stdio use m_GlobalSegMap, only : GlobalSegMap use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg use m_MCTWorld, only : ComponentToWorldRank use m_MCTWorld, only : ThisMCTWorld implicit noneINPUT PARAMETERS:
type(GlobalSegMap), intent(IN) :: outgoingGSMap integer, intent(IN) :: comp_id integer, intent(IN) :: TagBaseOUTPUT PARAMETERS:
integer, dimension(:), pointer :: reqHandle integer, optional, intent(OUT) :: statusREVISION HISTORY:
13Aug03 - J.W. Larson <[email protected]> - API and initial version. 05Mar04 - R. Jacob <[email protected]> - Send everything directly out of input GSMap. Don't use a SendBuffer.
This routine performs a blocking receive of a GlobalSegMap (the input argument outgoingGSMap) from the root processor on component comp_id. The input INTEGER argument TagBase is used to generate tags for the messages associated with this operation; there are six messages involved, so the user should avoid using tag values TagBase and TagBase + 5. The success (failure) of this operation is reported in the zero (non-zero) value of the optional INTEGER output variable status.
INTERFACE:
subroutine recv_(incomingGSMap, comp_id, TagBase, status)USES:
use m_mpif90 use m_die, only : MP_perr_die, die use m_stdio use m_GlobalSegMap, only : GlobalSegMap use m_GlobalSegMap, only : GlobalSegMap_init => init use m_MCTWorld, only : ComponentToWorldRank use m_MCTWorld, only : ThisMCTWorld implicit noneINPUT PARAMETERS:
integer, intent(IN) :: comp_id integer, intent(IN) :: TagBaseOUTPUT PARAMETERS:
type(GlobalSegMap), intent(OUT) :: incomingGSMap integer, optional, intent(OUT) :: statusREVISION HISTORY:
13Aug03 - J.W. Larson <[email protected]> - API and initial version. 25Aug03 - R.Jacob <[email protected]> - rename to recv_.
The routine bcast_() takes the input/output GlobalSegMap argument GSMap (on input valid only on the root process, on output valid on all processes) and broadcasts it to all processes on the communicator associated with the F90 handle comm. The success (failure) of this operation is returned as a zero (non-zero) value of the optional output INTEGER argument status.
INTERFACE:
subroutine bcast_(GSMap, root, comm, status)USES:
use m_mpif90 use m_die, only : MP_perr_die,die use m_stdio use m_GlobalSegMap, only : GlobalSegMap implicit noneINPUT PARAMETERS:
integer, intent(in) :: root integer, intent(in) :: commINPUT/OUTPUT PARAMETERS:
type(GlobalSegMap), intent(inout) :: GSMap ! Output GlobalSegMapOUTPUT PARAMETERS:
integer, optional, intent(out) :: status ! global vector sizeREVISION HISTORY:
17Oct01 - J.W. Larson <[email protected]> - Initial version. 11Aug03 - J.W. Larson <[email protected]> - Relocated from original location in m_GlobalSegMap.