next up previous contents
Next: 3 Global Segment Map Up: 1 Basic API's and Previous: 1 MCTWorld   Contents

Subsections

2 The Attribute Vector

2.1 Module m_AttrVect - Multi-field Storage (Source File: m_AttrVect.F90)

An attribute vector is a scheme for storing bundles of integer and real data vectors, indexed by the names of the fields stored in List format (see the mpeu module m_List for more information about the List datatype). The ordering of the fieldnames in the integer and real attribute List components (AttrVect%iList and AttrVect%rList, respectively) corresponds to the storage order of the attributes in their respective data buffers (the components AttrVect%iAttr(:,:) and AttrVect%rAttr(:,:), respectively). The organization of the fieldnames in List format, along with the direct mapping between List items and locations in the data buffer, allows the user to have random access to the field data. This approach also allows the user to set the number and the names of fields stored in an AttrVect at run-time.

The AttrVect stores field data in a pointwise fashion (that is, the data are grouped so that all the integer or real data associated with an individual point are adjacent to each other in memory. This amounts to the having the integer and real field data arrays in the AttrVect (the components AttrVect%iAttr(:,:) and AttrVect%rAttr(:,:), respectively) having the attribute index as the major (or fastest-varying) index. A prime example of this is observational data input to a data assimilation system. In the Model Coupling Toolkit, this datatype is the fundamental type for storing field data exchanged by component models, and forms a basis for other MCT datatypes that encapsulate time accumulation/averaging buffers (the Accumulator datatype defined in the module m_Accumulator), coordinate grid information (the GeneralGrid datatype defined in the module m_GeneralGrid), and sparse interpolation matrices (the SparseMatrix datatype defined in the module m_SparseMatrix).

The attribute vector is implemented in Fortran 90 using the AttrVect derived type. This module contains the definition of the AttrVect, and the numerous methods that service it. There are a number of initialization (creation) schemes, and a routine for zeroing out the elements of an AttrVect. There is a method to clean up allocated memory used by an AttrVect (destruction). There are numerous query methods that return: the number of datapoints (or length; the numbers of integer and real attributes; the data buffer index of a given real or integer attribute; and return the lists of real and integer attributes. There also exist methods for exporting a given attribute as a one-dimensional array and importing a given attribute from a one-dimensional array. There is a method for copying attributes from one AttrVect to another. There is also a method for cross-indexing the attributes in two AttrVect variables. In addition, there are methods that return those cross-indexed attributes along with some auxiliary data in a AVSharedIndicesOneType or AVSharedIndices structure. Finally, there are methods for sorting and permuting AttrVect entries using a MergeSort scheme keyed by the attributes of the AttrVect.


INTERFACE:

 
  module m_AttrVect
USES:
       use m_realkinds,only : SP,DP,FP          ! Real types definitions
 
       use m_List, only : List   ! Support for rList and iList components.
 
       implicit none
 
       private	! except
PUBLIC TYPES:
 
       public :: AttrVect        ! The class data structure
       public :: AVSharedIndicesOneType ! Data structure recording shared indices between
                                        ! two attribute vectors, for a single data type
                                        ! (e.g., shared real attributes)
       public :: AVSharedIndices ! Data structure recording shared indices between two
                                 ! attribute vectors, for all data types
 
     type AttrVect
 #ifdef SEQUENCE
       sequence
 #endif
       type(List) :: iList
       type(List) :: rList
       integer,dimension(:,:),pointer :: iAttr
       real(FP) ,dimension(:,:),pointer :: rAttr
     end type AttrVect
 
     type AVSharedIndicesOneType
        integer :: num_indices           ! number of shared items
        logical :: contiguous            ! true if index segments are contiguous in memory
        character*7 :: data_flag         ! data type flag (e.g., 'REAL' or 'INTEGER')
 
        ! arrays of indices to storage locations of shared attributes between the two
        ! attribute vectors: 
        integer, dimension(:), pointer :: aVindices1
        integer, dimension(:), pointer :: aVindices2
     end type AVSharedIndicesOneType
 
     type AVSharedIndices
        type(AVSharedIndicesOneType) :: shared_real     ! shared indices of type REAL
        type(AVSharedIndicesOneType) :: shared_integer  ! shared indices of type INTEGER
     end type AVSharedIndices
PUBLIC MEMBER FUNCTIONS:
 
       public :: init		! create a local vector
       public :: clean		! clean the local vector
       public :: zero		! zero the local vector
       public :: lsize		! size of the local vector
       public :: nIAttr		! number of integer attributes on local
       public :: nRAttr		! number of real attributes on local
       public :: indexIA		! index the integer attributes
       public :: indexRA		! index the real attributes
       public :: getIList        ! return list of integer attributes
       public :: getRList        ! return list of real attributes
       public :: exportIList     ! export INTEGER attibute List
       public :: exportRList     ! export REAL attibute List
       public :: exportIListToChar ! export INTEGER attibute List as Char
       public :: exportRListToChar ! export REAL attibute List as Char
       public :: appendIAttr     ! append INTEGER attribute List
       public :: appendRAttr     ! append REAL attribute List
       public :: exportIAttr     ! export INTEGER attribute to vector
       public :: exportRAttr     ! export REAL attribute to vector
       public :: importIAttr     ! import INTEGER attribute from vector
       public :: importRAttr     ! import REAL attribute from vector
       public :: Copy		! copy attributes from one Av to another
       public :: RCopy		! copy real attributes from one Av to another
       public :: ICopy		! copy integer attributes from one Av to another
       public :: Sort            ! sort entries, and return permutation
       public :: Permute         ! permute entries
       public :: Unpermute       ! Unpermute entries
       public :: SortPermute     ! sort and permute entries
       public :: SharedAttrIndexList  ! Cross-indices of shared
                                      ! attributes of two AttrVects
       public :: SharedIndices        ! Given two AttrVects, create an AVSharedIndices structure
       public :: SharedIndicesOneType ! Given two AttrVects, create an
                                      ! AVSharedIndicesOneType structure for a single type
       public :: cleanSharedIndices   ! clean a AVSharedIndices structure
       public :: cleanSharedIndicesOneType  ! clean a AVSharedIndicesOneType structure
 
 
     interface init   ; module procedure	&
        init_,  &
        initv_, &
        initl_
     end interface
     interface clean  ; module procedure clean_  ; end interface
     interface zero  ; module procedure zero_  ; end interface
     interface lsize  ; module procedure lsize_  ; end interface
     interface nIAttr ; module procedure nIAttr_ ; end interface
     interface nRAttr ; module procedure nRAttr_ ; end interface
     interface indexIA; module procedure indexIA_; end interface
     interface indexRA; module procedure indexRA_; end interface
     interface getIList; module procedure getIList_; end interface
     interface getRList; module procedure getRList_; end interface
     interface exportIList; module procedure exportIList_; end interface
     interface exportRList; module procedure exportRList_; end interface
     interface exportIListToChar
        module procedure exportIListToChar_
     end interface
     interface exportRListToChar
        module procedure exportRListToChar_
     end interface
     interface appendIAttr  ; module procedure appendIAttr_  ; end interface
     interface appendRAttr  ; module procedure appendRAttr_  ; end interface
     interface exportIAttr; module procedure exportIAttr_; end interface
     interface exportRAttr; module procedure &
        exportRAttrSP_, &
        exportRAttrDP_
     end interface
     interface importIAttr; module procedure importIAttr_; end interface
     interface importRAttr; module procedure &
          importRAttrSP_, &
          importRAttrDP_
     end interface
     interface Copy    ; module procedure Copy_    ; end interface
     interface RCopy   ; module procedure  &
          RCopy_, &
          RCopyL_
     end interface
     interface ICopy   ; module procedure  &
           ICopy_, &
           ICopyL_
     end interface
     interface Sort    ; module procedure Sort_    ; end interface
     interface Permute ; module procedure Permute_ ; end interface
     interface Unpermute ; module procedure Unpermute_ ; end interface
     interface SortPermute ; module procedure SortPermute_ ; end interface
     interface SharedAttrIndexList ; module procedure &
         aVaVSharedAttrIndexList_ 
     end interface
     interface SharedIndices ; module procedure SharedIndices_ ; end interface
     interface SharedIndicesOneType ; module procedure SharedIndicesOneType_ ; end interface
     interface cleanSharedIndices ; module procedure cleanSharedIndices_ ; end interface
     interface cleanSharedIndicesOneType ; module procedure cleanSharedIndicesOneType_ ; end interface
REVISION HISTORY:
   10Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   10Oct00 - J.W. Larson <[email protected]> - made getIList
             and getRList functions public and added appropriate
             interface definitions
   20Oct00 - J.W. Larson <[email protected]> - added Sort, 
             Permute, and SortPermute functions.
   09May01 - J.W. Larson <[email protected]> - added initl_().
   19Oct01 - J.W. Larson <[email protected]> - added routines 
             exportIattr(), exportRAttr(), importIAttr(), 
             and importRAttr().  Also cleaned up module and 
             routine prologues.
   13Dec01 - J.W. Larson <[email protected]> - made importIAttr()
             and importRAttr() public (bug fix).
   14Dec01 - J.W. Larson <[email protected]> - added exportIList()
             and exportRList().
   14Feb02 - J.W. Larson <[email protected]> - added CHARCTER
             functions exportIListToChar() and exportRListToChar()
   26Feb02 - J.W. Larson <[email protected]> - corrected of usage
             of m_die routines throughout this module.
   16Apr02 - J.W. Larson <[email protected]> - added the method
             LocalReduce(), and the public data members AttrVectSUM,
             AttrVectMIN, and AttrVectMAX.
   7May02 - J.W. Larson <[email protected]> - Refactoring.  Moved
            LocalReduce() and the public data members AttrVectSUM,
             AttrVectMIN, and AttrVectMAX to a new module named
             m_AttrVectReduce.
   12Jun02 - R.L. Jacob <[email protected]> - add Copy function
   13Jun02 - R.L. Jacob <[email protected]> - move aVavSharedAttrIndexList
             to this module from old m_SharedAttrIndicies
   28Apr11 - W.J. Sacks <[email protected]> - added AVSharedIndices and
             AVSharedIndicesOneType derived types, and associated
             subroutines
   10Apr12 - W.J. Sacks <[email protected]> - modified AVSharedIndices code
             to be Fortran-90 compliant

2.1.1 init_ - Initialize an AttrVect Given Attribute Lists and Length

This routine creates an AttrVect (the output argument aV) using the optional input CHARACTER arguments iList, and rList to define its integer and real attributes, respectively. The optional input INTEGER argument lsize defines the number of points for which we are storing attributes, or the length of aV. The expected form for the arguments iList and rList are colon-delimited strings where each substring defines an attribute. Suppose we wish to store N observations that have the real attributes 'latitude', 'longitude', pressure, 'u-wind', and 'v-wind'. Suppose we also wish to store the integer attributes 'hour', 'day', 'month', 'year', and 'data source'. This can be accomplished by invoking init_() as follows:

 
   call init_(aV, 'hour:day:month:year:data source', &
              'latitude:longitude:pressure:u-wind:v-wind', N)
The resulting AttrVect aV will have five integer attributes, five real attributes, and length N.


INTERFACE:

 
  subroutine init_(aV, iList, rList, lsize)
USES:
       use m_List, only : List
       use m_List, only : init,nitem
       use m_List, only : List_nullify => nullify
       use m_mall
       use m_die
 
       implicit none
INPUT PARAMETERS:
       character(len=*), optional, intent(in)  :: iList
       character(len=*), optional, intent(in)  :: rList
       integer,          optional, intent(in)  :: lsize
OUTPUT PARAMETERS:
       type(AttrVect),             intent(out) :: aV
REVISION HISTORY:
   09Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   09Oct01 - J.W. Larson <[email protected]> - added feature to
             nullify all pointers before usage.  This was done to
             accomodate behavior of the f90 ASSOCIATED intrinsic 
             function on the AIX platform.
   07Dec01 - E.T. Ong <[email protected]> - added support for 
             intialization with blank character strings for iList
             and rList

2.1.2 initv_ - Initialize One AttrVect from Another

This routine takes an input AttrVect argument bV, and uses its attribute list information to create an output AttrVect variable aV. The length of aV is defined by the input INTEGER argument lsize.


INTERFACE:

 
  subroutine initv_(aV, bV, lsize)
USES:
       use m_String, only : String,char
       use m_String, only : String_clean => clean    
       use m_List,   only : get
       use m_List,   only : List_nullify => nullify
       use m_die
       use m_stdio
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),intent(in)  :: bV
       integer,       intent(in)  :: lsize
OUTPUT PARAMETERS:
       type(AttrVect),intent(out) :: aV
REVISION HISTORY:
   22Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   17May01 - R. Jacob <[email protected]> - add a check to see if
             input argument has been defined.  SGI will dump
             core if its not.
   10Oct01 - J. Larson <[email protected]> - Nullify all pointers
             in ouput AttrVect aV before initializing aV.
   19Sep08 - J. Wolfe <[email protected]> - plug memory leak from not deallocating
             strings.

2.1.3 initl_ - Initialize an AttrVect Using the List Type

This routine initializes an AttrVect directly from input List data type arguments iList and rList (see the module m_List in mpeu for further details), and an input length lsize. The resulting AttrVect is returned in the argument aV.

N.B.: If the user supplies an empty list for the arguments iList (rList), then aV will be created only with REAL (INTEGER) attributes. If both arguments iList and rList are empty, the routine will terminate execution and report an error.

N.B.: The outcome of this routine, aV represents allocated memory. When this AttrVect is no longer needed, it must be deallocated by invoking the routine AttrVect_clean(). Failure to do so will spawn a memory leak.


INTERFACE:

 
  subroutine initl_(aV, iList, rList, lsize)
USES:
       use m_die
       use m_stdio
 
       use m_String, only : String
       use m_String, only : String_clean => clean
       use m_String, only : String_toChar => toChar
 
       use m_List, only : List
       use m_List, only : List_nitem => nitem
       use m_List, only : List_exportToChar => exportToChar
 
       implicit none
INPUT PARAMETERS:
       type(List),  intent(in)  :: iList
       type(List),  intent(in)  :: rList
       integer,     intent(in)  :: lsize
OUTPUT PARAMETERS:
       type(AttrVect), intent(out) :: aV
REVISION HISTORY:
   09May98 - J.W. Larson <[email protected]> - initial version.
   08Aug01 - E.T. Ong <[email protected]> - change list assignment(=)
             to list copy to avoid compiler errors with pgf90.
   10Oct01 - J. Larson <[email protected]> - Nullify all pointers
             in ouput AttrVect aV before initializing aV.  Also, 
             greater caution taken regarding validity of input 
             arguments iList and rList.
   15May08 - J. Larson <[email protected]> - Simplify to use
             the init_ routine.  Better argument checking.

2.1.4 clean_ - Deallocate Allocated Memory Structures of an AttrVect

This routine deallocates the allocated memory structures of the input/output AttrVect argument aV. This amounts to cleaning the List structures aV%iList and av%rList, and deallocating the arrays aV%iAttr(:,:) and aV%rAttr(:,:). The success (failure) of this operation is signified by a zero (non-zero) value of the optional INTEGER output argument stat. If clean_() is invoked without supplying stat, and any of the deallocation operations fail, the routine will terminate with an error message.


INTERFACE:

 
  subroutine clean_(aV, stat)
USES:
       use m_mall
       use m_stdio
       use m_die
       use m_List, only : List_clean => clean
 
       implicit none
INPUT/OUTPUT PARAMETERS:
       type(AttrVect),    intent(inout) :: aV
OUTPUT PARAMETERS:
       integer, optional, intent(out)   :: stat
REVISION HISTORY:
   09Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   10Oct01 - J. Larson <[email protected]> - various fixes to 
             prevent deallocation of UNASSOCIATED pointers.
   01Mar01 - E.T. Ong <[email protected]> - removed dies to prevent
             crashes when cleaning uninitialized attrvects. Added
             optional stat argument.

2.1.5 lsize_ - Length of an AttrVect

This function returns the number of elements, or length of the input AttrVect argument aV. This function examines the length of the second dimension of the arrays aV%iAttr(:,:) and aV%rAttr(:,:). If neither aV%iAttr(:,:) nor aV%rAttr(:,:) are associated, then ${\tt lsize\_(aV)} = 0$. If aV%iAttr(:,:) is associated, but aV%rAttr(:,:) is not, then ${\tt lsize\_(aV)} = {\tt size(aV\%iAttr,2)}$. If aV%iAttr(:,:) is not associated, but aV%rAttr(:,:) is, then ${\tt lsize\_(aV)} = {\tt size(aV\%rAttr,2)}$. If both aV%iAttr(:,:) and aV%rAttr(:,:) are associated, the function lsize_() will do one of two things: If ${\tt size(aV\%iAttr,2)} = {\tt size(aV\%rAttr,2)}$, this equal value will be returned. If ${\tt size(aV\%iAttr,2)} \neq
{\tt size(aV\%rAttr,2)}$, termination with an error message will occur.


INTERFACE:

 
  integer function lsize_(aV)
USES:
 
      use m_List,  only : List
      use m_List,  only : List_allocated => allocated
 
      use m_stdio, only : stderr
      use m_die
  
      implicit none
INPUT PARAMETERS:
       type(AttrVect), intent(in) :: aV
REVISION HISTORY:
   09Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   10Oct01 - J. Larson <[email protected]> - made code more robust
             to handle cases where the length of either aV%iAttr or
             aV%rAttr is zero, but the other is positive.\end{verbatim}
 
%/////////////////////////////////////////////////////////////
 
\mbox{}\hrulefill\ 
 

  \subsubsection{zero\_ - Set AttrVect Field Data to Zero}


  
   This routine sets all of the point values of the integer and real 
   attributes of an the input/output {\tt AttrVect} argument {\tt aV}
   to zero.  The default action is to set the values of all the real and 
   integer attributes to zero.  The user may prevent the zeroing of the 
   real (integer) attributes invoking {\tt zero\_()} with the optional 
   {\tt LOGICAL} argument {\tt zeroReals} ({\tt zeroInts}) set with value 
   {\tt .FALSE.}
  
\bigskip{\sf INTERFACE:}
\begin{verbatim} 
  subroutine zero_(aV, zeroReals, zeroInts)
USES:
 
 
      use m_die,only	: die
      use m_stdio,only	: stderr
 
      use m_List, only : List
      use m_List, only : List_allocated => allocated
  
      implicit none
INPUT PARAMETERS:
 
      logical, optional, intent(IN)    :: zeroReals
      logical, optional, intent(IN)    :: zeroInts
INPUT/OUTPUT PARAMETERS:
      type(AttrVect),    intent(INOUT) :: aV
REVISION HISTORY:
   17May01 - R. Jacob <[email protected]> - initial prototype/code
   15Oct01 - J. Larson <[email protected]> - switched loop order
             for cache optimization.
   03Dec01 - E.T. Ong <[email protected]> - eliminated looping method of
             of zeroing. "Compiler assignment" of attrvect performs faster
             on the IBM SP with mpxlf90 compiler.
   05Jan10 - R. Jacob <[email protected]> - zeroing an uninitialized aV is no
             longer a fatal error.

2.1.6 nIAttr_ - Return the Number of Integer Attributes

This integer function returns the number of integer attributes present in the input AttrVect argument aV.


INTERFACE:

 
  integer function nIAttr_(aV)
USES:
       use m_List, only : nitem
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),intent(in) :: aV
REVISION HISTORY:
   22Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   10Oct01 - J. Larson <[email protected]> - made code more robust
             by checking status of pointers in aV%iList\end{verbatim}
 
%/////////////////////////////////////////////////////////////
 
\mbox{}\hrulefill\ 
 

  \subsubsection{nRAttr\_ - Return the Number of Real Attributes}


  
   This integer function returns the number of real attributes 
   present in the input {\tt AttrVect} argument {\tt aV}.
 
\bigskip{\sf INTERFACE:}
\begin{verbatim} 
  integer function nRAttr_(aV)
USES:
       use m_List, only : nitem
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),intent(in) :: aV
REVISION HISTORY:
   22Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   10Oct01 - J. Larson <[email protected]> - made code more robust
             by checking status of pointers in aV%iList\end{verbatim}
 
%/////////////////////////////////////////////////////////////
 
\mbox{}\hrulefill\ 
 

  \subsubsection{getIList\_ - Retrieve the Name of a Numbered Integer Attribute}


  
   This routine returns the name of the {\tt ith} integer attribute of 
   the input {\tt AttrVect} argument {\tt aVect}.  The name is returned 
   in the output {\tt String} argument {\tt item} (see the mpeu module 
   {\tt m\_String} for more information regarding the {\tt String} type).
  
\bigskip{\sf INTERFACE:}
\begin{verbatim} 
  subroutine getIList_(item, ith, aVect)
USES:
       use m_String, only : String
       use m_List,   only : get
 
       implicit none
INPUT PARAMETERS:
       integer,     intent(in)  :: ith
       type(AttrVect),intent(in) :: aVect
OUTPUT PARAMETERS:
       type(String),intent(out) :: item
REVISION HISTORY:
   24Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code

2.1.7 getRList_ - Retrieve the Name of a Numbered Real Attribute

This routine returns the name of the ith real attribute of the input AttrVect argument aVect. The name is returned in the output String argument item (see the mpeu module m_String for more information regarding the String type).


INTERFACE:

 
  subroutine getRList_(item, ith, aVect)
USES:
       use m_String, only : String
       use m_List,   only : get
 
       implicit none
INPUT PARAMETERS:
       integer,        intent(in)  :: ith
       type(AttrVect), intent(in)  :: aVect
OUTPUT PARAMETERS:
       type(String),   intent(out) :: item
REVISION HISTORY:
   24Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code

2.1.8 indexIA_ - Index an Integer Attribute

This function returns an INTEGER, corresponding to the location of an integer attribute within the input AttrVect argument aV. For example, suppose aV has the following attributes 'month', 'day', and 'year'. The array of integer values for the attribute 'day' is stored in aV%iAttr(indexIA_(aV,'day'),:). If indexIA_() is unable to match item to any of the integer attributes in aV, the resulting value is zero which is equivalent to an error. The optional input CHARACTER arguments perrWith and dieWith control how such errors are handled.

  1. if neither perrWith nor dieWith are present, indexIA_() terminates execution with an internally generated error message;
  2. if perrWith is present, but dieWith is not, an error message is written to stderr incorporating user-supplied traceback information stored in the argument perrWith;
  3. if perrWith is present, but dieWith is not, and perrWith is equal to ``quiet'', no error message is written.
  4. if dieWith is present, execution terminates with an error message written to stderr that incorporates user-supplied traceback information stored in the argument dieWith; and
  5. if both perrWith and dieWith are present, execution terminates with an error message using dieWith, and the argument perrWith is ignored.


INTERFACE:

 
  integer function indexIA_(aV, item, perrWith, dieWith)
USES:
       use m_die,  only : die
       use m_stdio,only : stderr
 
       use m_String, only : String
       use m_String, only : String_init => init
       use m_String, only : String_clean => clean
       use m_String, only : String_ToChar => ToChar
 
       use m_List, only : index
 
       use m_TraceBack, only : GenTraceBackString
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),             intent(in) :: aV
       character(len=*),           intent(in) :: item
       character(len=*), optional, intent(in) :: perrWith
       character(len=*), optional, intent(in) :: dieWith
REVISION HISTORY:
   27Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
    2Aug02 - J. Larson - Solidified error handling using perrWith/dieWith
    1Jan05 - R. Jacob - add quiet option for error handling

2.1.9 indexRA_ - Index a Real Attribute

This function returns an INTEGER, corresponding to the location of a real attribute within the input AttrVect argument aV. For example, suppose aV has the following attributes 'latitude', 'longitude', and 'pressure'. The array of real values for the attribute 'longitude' is stored in aV%iAttr(indexRA_(aV,'longitude'),:). If indexRA_() is unable to match item to any of the real attributes in aV, the resulting value is zero which is equivalent to an error. The optional input CHARACTER arguments perrWith and dieWith control how such errors are handled.

  1. if neither perrWith nor dieWith are present, indexRA_() terminates execution with an internally generated error message;
  2. if perrWith is present, but dieWith is not, an error message is written to stderr incorporating user-supplied traceback information stored in the argument perrWith;
  3. if perrWith is present, but dieWith is not, and perrWith is equal to ``quiet'', no error message is written.
  4. if dieWith is present, execution terminates with an error message written to stderr that incorporates user-supplied traceback information stored in the argument dieWith; and
  5. if both perrWith and dieWith are present, execution terminates with an error message using dieWith, and the argument perrWith is ignored.


INTERFACE:

 
  integer function indexRA_(aV, item, perrWith, dieWith)
USES:
       use m_die,  only : die
       use m_stdio,only : stderr
 
       use m_List, only : index
 
       use m_String, only : String
       use m_String, only : String_init => init
       use m_String, only : String_clean => clean
       use m_String, only : String_ToChar => ToChar
 
       use m_TraceBack, only : GenTraceBackString
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),             intent(in) :: aV
       character(len=*),           intent(in) :: item
       character(len=*), optional, intent(in) :: perrWith
       character(len=*), optional, intent(in) :: dieWith
REVISION HISTORY:
   27Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
    2Aug02 - J. Larson - Solidified error handling using perrWith/dieWith
   18Jan05 - R. Jacob - add quiet option for error handling

2.1.10 appendIAttr_ - Append one or more attributes onto the INTEGER part of an AttrVect.

This routine takes an input AttrVect argument aV, and an input character string rList and Appends rList to the INTEGER part of aV. The success (failure) of this operation is signified by a zero (nonzero) value for the optional INTEGER output argument status.


INTERFACE:

 
  subroutine appendIAttr_(aV, iList, status)
USES:
       use m_List,   only : List_init => init
       use m_List,   only : List_append => append
       use m_List,   only : List_clean => clean
       use m_List,   only : List_nullify => nullify
       use m_List,   only : List_allocated => allocated
       use m_List,   only : List_copy => copy
       use m_List,   only : List
       use m_die
       use m_stdio
 
       implicit none
INPUT/OUTPUT PARAMETERS:
       type(AttrVect),intent(inout)  :: aV
INPUT PARAMETERS:
       character(len=*), intent(in)  :: iList
OUTPUT PARAMETERS:
       integer,optional,intent(out)  :: status
REVISION HISTORY:
   08Jul03 - R. Jacob <[email protected]> - initial version

2.1.11 appendRAttr_ - Append one or more attributes onto the REAL part of an AttrVect.

This routine takes an input AttrVect argument aV, and an input character string rList and Appends rList to the REAL part of aV. The success (failure) of this operation is signified by a zero (nonzero) value for the optional INTEGER output argument status.


INTERFACE:

 
  subroutine appendRAttr_(aV, rList, status)
USES:
       use m_List,   only : List_init => init
       use m_List,   only : List_append => append
       use m_List,   only : List_clean => clean
       use m_List,   only : List_nullify => nullify
       use m_List,   only : List_allocated => allocated
       use m_List,   only : List_copy => copy
       use m_List,   only : List
       use m_die
       use m_stdio
 
       implicit none
INPUT/OUTPUT PARAMETERS:
       type(AttrVect),intent(inout)  :: aV
INPUT PARAMETERS:
       character(len=*), intent(in)  :: rList
OUTPUT PARAMETERS:
       integer,optional,intent(out)  :: status
REVISION HISTORY:
   04Jun03 - R. Jacob <[email protected]> - initial version

2.1.12 exportIList_ - Return INTEGER Attribute List

This routine extracts from the input AttrVect argument aV the integer attribute list, and returns it as the List output argument outIList. The success (failure) of this operation is signified by a zero (nonzero) value for the optional INTEGER output argument status.

N.B.: This routine returns an allocated List data structure (outIList). The user is responsible for deallocating this structure by invoking List_clean() (see the module m_List for details) once it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine exportIList_(aV, outIList, status)
USES:
       use m_die ,  only : die
       use m_stdio, only : stderr
 
       use m_List,  only : List
       use m_List,  only : List_allocated => allocated
       use m_List,  only : List_copy => copy
       use m_List,  only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in)  :: aV
OUTPUT PARAMETERS:
 
       type(List),                 intent(out) :: outIList
       integer,          optional, intent(out) :: status
REVISION HISTORY:
   14Dec01 - J.W. Larson <[email protected]> - initial prototype.

2.1.13 exportRList_ - Return REAL attribute List

This routine extracts from the input AttrVect argument aV the real attribute list, and returns it as the List output argument outRList. The success (failure) of this operation is signified by a zero (nonzero) value for the optional INTEGER output argument status.

N.B.: This routine returns an allocated List data structure (outRList). The user is responsible for deallocating this structure by invoking List_clean() (see the module m_List for details) once it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine exportRList_(aV, outRList, status)
USES:
       use m_die ,  only : die
       use m_stdio, only : stderr
 
       use m_List,  only : List
       use m_List,  only : List_allocated => allocated
       use m_List,  only : List_copy => copy
       use m_List,  only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),           intent(in)  :: aV
OUTPUT PARAMETERS:
 
       type(List),               intent(out) :: outRList
       integer,        optional, intent(out) :: status
REVISION HISTORY:
   14Dec01 - J.W. Larson <[email protected]> - initial prototype.

2.1.14 exportIListToChar_ - Return AttrVect%iList as CHARACTER

This routine extracts from the input AttrVect argument aV the integer attribute list (see the mpeu module m_List for more information regarding the List type), and returns it as a CHARACTER suitable for printing. An example of its usage is

             write(stdout,'(1a)') exportIListToChar_(aV)
which writes the contents of aV%iList%bf to the Fortran device stdout.


INTERFACE:

 
  function exportIListToChar_(aV)
USES:
       use m_die ,  only : die
       use m_stdio, only : stderr
 
       use m_List,  only : List
       use m_List,  only : List_allocated => allocated
       use m_List,  only : List_copy => copy
       use m_List,  only : List_exportToChar => exportToChar
       use m_List,  only : List_clean => clean
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),       intent(in) :: aV
OUTPUT PARAMETERS:
 
       character(len=size(aV%iList%bf,1)) :: exportIListToChar_
 
REVISION HISTORY:
   13Feb02 - J.W. Larson <[email protected]> - initial prototype.
   05Jun03 - R. Jacob <[email protected]> - return a blank instead of dying
             to avoid I/O errors when this function is used in a write statement.

2.1.15 exportRListToChar_ - Return AttrVect%rList as CHARACTER

This routine extracts from the input AttrVect argument aV the real attribute list (see the mpeu module m_List for more information regarding the List type), and returns it as a CHARACTER suitable for printing. An example of its usage is

             write(stdout,'(1a)') exportRListToChar_(aV)
which writes the contents of aV%rList%bf to the Fortran device stdout.


INTERFACE:

 
  function exportRListToChar_(aV)
USES:
       use m_die ,  only : die
       use m_stdio, only : stderr
 
       use m_List,  only : List
       use m_List,  only : List_allocated => allocated
       use m_List,  only : List_copy => copy
       use m_List,  only : List_exportToChar => exportToChar
       use m_List,  only : List_clean => clean
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),       intent(in) :: aV
OUTPUT PARAMETERS:
 
       character(len=size(aV%rList%bf,1)) :: exportRListToChar_
 
REVISION HISTORY:
   13Feb02 - J.W. Larson <[email protected]> - initial prototype.
   05Jun03 - R. Jacob <[email protected]> - return a blank instead of dying
             to avoid I/O errors when this function is used in a write statement.

2.1.16 exportIAttr_ - Return INTEGER Attribute as a Vector

This routine extracts from the input AttrVect argument aV the integer attribute corresponding to the tag defined in the input CHARACTER argument AttrTag, and returns it in the INTEGER output array outVect, and its length in the output INTEGER argument lsize. The optional input CHARACTER arguments perrWith and dieWith control how errors are handled.

  1. if neither perrWith nor dieWith are present, exportIAttr_() terminates execution with an internally generated error message;
  2. if perrWith is present, but dieWith is not, an error message is written to stderr incorporating user-supplied traceback information stored in the argument perrWith;
  3. if dieWith is present, execution terminates with an error message written to stderr that incorporates user-supplied traceback information stored in the argument dieWith; and
  4. if both perrWith and dieWith are present, execution terminates with an error message using dieWith, and the argument perrWith is ignored.

N.B.: This routine will fail if the AttrTag is not in the AttrVect List component aV%iList.

N.B.: The flexibility of this routine regarding the pointer association status of the output argument outVect means the user must invoke this routine with care. If the user wishes this routine to fill a pre-allocated array, then obviously this array must be allocated prior to calling this routine. If the user wishes that the routine create the output argument array outVect, then the user must ensure this pointer is not allocated (i.e. the user must nullify this pointer) before this routine is invoked.

N.B.: If the user has relied on this routine to allocate memory associated with the pointer outVect, then the user is responsible for deallocating this array once it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine exportIAttr_(aV, AttrTag, outVect, lsize, perrWith, dieWith)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       use m_String, only : String
       use m_String, only : String_init => init
       use m_String, only : String_clean => clean
       use m_String, only : String_ToChar => ToChar
 
       use m_TraceBack, only : GenTraceBackString
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in) :: aV
       character(len=*),           intent(in) :: AttrTag
       character(len=*), optional, intent(in) :: perrWith
       character(len=*), optional, intent(in) :: dieWith
OUTPUT PARAMETERS:
 
       integer,      dimension(:), pointer     :: outVect
       integer,          optional, intent(out) :: lsize
REVISION HISTORY:
   19Oct01 - J.W. Larson <[email protected]> - initial (slow) 
             prototype.
    6May02 - J.W. Larson <[email protected]> - added capability 
             to work with pre-allocated outVect.

2.1.17 exportRAttrSP_ - Return REAL Attribute as a Pointer to Array

This routine extracts from the input AttrVect argument aV the real attribute corresponding to the tag defined in the input CHARACTER argument AttrTag, and returns it in the REAL output array outVect, and its length in the output INTEGER argument lsize. The optional input CHARACTER arguments perrWith and dieWith control how errors are handled.

  1. if neither perrWith nor dieWith are present, exportRAttr_() terminates execution with an internally generated error message;
  2. if perrWith is present, but dieWith is not, an error message is written to stderr incorporating user-supplied traceback information stored in the argument perrWith;
  3. if dieWith is present, execution terminates with an error message written to stderr that incorporates user-supplied traceback information stored in the argument dieWith; and
  4. if both perrWith and dieWith are present, execution terminates with an error message using dieWith, and the argument perrWith is ignored.

N.B.: This routine will fail if the AttrTag is not in the AttrVect List component aV%iList.

N.B.: The flexibility of this routine regarding the pointer association status of the output argument outVect means the user must invoke this routine with care. If the user wishes this routine to fill a pre-allocated array, then obviously this array must be allocated prior to calling this routine. If the user wishes that the routine create the output argument array outVect, then the user must ensure this pointer is not allocated (i.e. the user must nullify this pointer) before this routine is invoked.

N.B.: If the user has relied on this routine to allocate memory associated with the pointer outVect, then the user is responsible for deallocating this array once it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine exportRAttrSP_(aV, AttrTag, outVect, lsize, perrWith, dieWith)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
 
       use m_String, only : String
       use m_String, only : String_init => init
       use m_String, only : String_clean => clean
       use m_String, only : String_ToChar => ToChar
 
       use m_TraceBack, only : GenTraceBackString
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in) :: aV
       character(len=*),           intent(in) :: AttrTag
       character(len=*), optional, intent(in) :: perrWith
       character(len=*), optional, intent(in) :: dieWith
OUTPUT PARAMETERS:
 
       real(SP),        dimension(:),  pointer     :: outVect
       integer,          optional,     intent(out) :: lsize
REVISION HISTORY:
   19Oct01 - J.W. Larson <[email protected]> - initial (slow) 
             prototype.
    6May02 - J.W. Larson <[email protected]> - added capability 
             to work with pre-allocated outVect.

2.1.18 importIAttr_ - Import INTEGER Vector as an Attribute

This routine imports into the input/output AttrVect argument aV the integer attribute corresponding to the tag defined in the input CHARACTER argument AttrTag. The data to be imported is provided in the INTEGER input array inVect, and the number of entries to be imported in the optional input INTEGER argument lsize.

N.B.: This routine will fail if the AttrTag is not in the AttrVect List component aV%iList.


INTERFACE:

 
  subroutine importIAttr_(aV, AttrTag, inVect, lsize)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       implicit none
INPUT PARAMETERS:
 
       character(len=*),       intent(in)    :: AttrTag
       integer,  dimension(:), pointer       :: inVect
       integer,  optional,     intent(in)    :: lsize
INPUT/OUTPUT PARAMETERS:
 
       type(AttrVect),         intent(inout) :: aV
REVISION HISTORY:
   19Oct01 - J.W. Larson <[email protected]> - initial (slow) 
             prototype.

2.1.19 importRAttrSP_ - Import REAL Vector as an Attribute

This routine imports into the input/output AttrVect argument aV the real attribute corresponding to the tag defined in the input CHARACTER argument AttrTag. The data to be imported is provided in the REAL input array inVect, and its length in the optional input INTEGER argument lsize.

N.B.: This routine will fail if the AttrTag is not in the AttrVect List component aV%rList.


INTERFACE:

 
  subroutine importRAttrSP_(aV, AttrTag, inVect, lsize)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       implicit none
INPUT PARAMETERS:
 
       character(len=*),   intent(in)    :: AttrTag
       real(SP), dimension(:), pointer   :: inVect
       integer, optional,  intent(in)    :: lsize
INPUT/OUTPUT PARAMETERS:
 
       type(AttrVect),     intent(inout) :: aV
REVISION HISTORY:
   19Oct01 - J.W. Larson <[email protected]> - initial (slow) 
             prototype.

2.1.20 RCopy_ - Copy Real Attributes from One AttrVect to Another

This routine copies from input argment aVin into the output AttrVect argument aVout the shared real attributes.

If the optional argument Vector is present and true, the vector architecture-friendly portions of this routine will be invoked.

If the optional argument sharedIndices is present, it should be the result of the call SharedIndicesOneType_(aVin, aVout, 'REAL', sharedIndices). Providing this argument speeds up this routine substantially. For example, you can compute a sharedIndices structure once for a given pair of AttrVects, then use that same structure for all copies between those two AttrVects (although note that a different sharedIndices variable would be needed if aVin and aVout were reversed).

N.B.: This routine will fail if the aVout is not initialized.


INTERFACE:

 
  subroutine RCopy_(aVin, aVout, vector, sharedIndices)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       use m_List,          only : GetSharedListIndices
       use m_List,          only : GetIndices => get_indices
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in)    :: aVin
       logical, optional,          intent(in)    :: vector 
       type(AVSharedIndicesOneType), optional, intent(in) :: sharedIndices
OUTPUT PARAMETERS:
 
       type(AttrVect),             intent(inout) :: aVout
REVISION HISTORY:
   18Aug06 - R. Jacob <[email protected]> - initial version.
   28Apr11 - W.J. Sacks <[email protected]> - added sharedIndices argument

2.1.21 RCopyL_ - Copy Specific Real Attributes from One AttrVect to Another

This routine copies from input argment aVin into the output AttrVect argument aVout the real attributes specified in input CHARACTER argument rList. The attributes can be listed in any order.

If any attributes in aVout have different names but represent the the same quantity and should still be copied, you must provide a translation argument TrList. The translation arguments should be identical in length to the rList but with the correct aVout name subsititued at the appropriate place.

If the optional argument Vector is present and true, the vector architecture-friendly portions of this routine will be invoked.

N.B.: This routine will fail if the aVout is not initialized or if any of the specified attributes are not present in either aVout or aVin.


INTERFACE:

 
  subroutine RCopyL_(aVin, aVout, rList, TrList,  vector)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       use m_List,          only : GetSharedListIndices
       use m_List,          only : GetIndices => get_indices
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in)    :: aVin
       character(len=*),           intent(in)    :: rList
       character(len=*), optional, intent(in)    :: TrList
       logical, optional,          intent(in)    :: vector
OUTPUT PARAMETERS:
 
       type(AttrVect),             intent(inout) :: aVout
REVISION HISTORY:
   16Aug06 - R. Jacob <[email protected]> - initial version from breakup
             of Copy_.

2.1.22 ICopy_ - Copy Integer Attributes from One AttrVect to Another

This routine copies from input argment aVin into the output AttrVect argument aVout the shared integer attributes.

If the optional argument Vector is present and true, the vector architecture-friendly portions of this routine will be invoked.

If the optional argument sharedIndices is present, it should be the result of the call SharedIndicesOneType_(aVin, aVout, 'INTEGER', sharedIndices). Providing this argument speeds up this routine substantially. For example, you can compute a sharedIndices structure once for a given pair of AttrVects, then use that same structure for all copies between those two AttrVects (although note that a different sharedIndices variable would be needed if aVin and aVout were reversed).

N.B.: This routine will fail if the aVout is not initialized.


INTERFACE:

 
  subroutine ICopy_(aVin, aVout, vector, sharedIndices)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       use m_List,          only : GetSharedListIndices
       use m_List,          only : GetIndices => get_indices
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in)    :: aVin
       logical, optional,          intent(in)    :: vector 
       type(AVSharedIndicesOneType), optional, intent(in) :: sharedIndices
OUTPUT PARAMETERS:
 
       type(AttrVect),             intent(inout) :: aVout
REVISION HISTORY:
   16Aug06 - R. Jacob <[email protected]> - initial version.
   28Apr11 - W.J. Sacks <[email protected]> - added sharedIndices argument

2.1.23 ICopyL_ - Copy Specific Integer Attributes from One AttrVect to Another

This routine copies from input argment aVin into the output AttrVect argument aVout the integer attributes specified in input CHARACTER argument iList. The attributes can be listed in any order.

If any attributes in aVout have different names but represent the the same quantity and should still be copied, you must provide a translation argument TiList. The translation arguments should be identical in length to the iList but with the correct aVout name subsititued at the appropriate place.

If the optional argument Vector is present and true, the vector architecture-friendly portions of this routine will be invoked.

N.B.: This routine will fail if the aVout is not initialized or if any of the specified attributes are not present in either aVout or aVin.


INTERFACE:

 
  subroutine ICopyL_(aVin, aVout, iList, TiList, vector)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       use m_List,          only : GetIndices => get_indices
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in)    :: aVin
       character(len=*)          , intent(in)    :: iList
       character(len=*), optional, intent(in)    :: TiList
       logical, optional,          intent(in)    :: vector
OUTPUT PARAMETERS:
 
       type(AttrVect),             intent(inout) :: aVout
REVISION HISTORY:
   16Aug06 - R. Jacob <[email protected]> - initial version from breakup
             of Copy_.

2.1.24 Copy_ - Copy Real and Integer Attributes from One AttrVect to Another

This routine copies from input argment aVin into the output AttrVect argument aVout the real and integer attributes specified in input CHARACTER argument iList and rList. The attributes can be listed in any order. If neither iList nor rList are provided, all attributes shared between aVin and aVout will be copied.

If any attributes in aVout have different names but represent the the same quantity and should still be copied, you must provide a translation argument TrList and/or TiList. The translation arguments should be identical to the rList or iList but with the correct aVout name subsititued at the appropriate place.

This routines combines the functions of RCopy_, RCopyL_, ICopy_ and ICopyL_. If you know you only want to copy real attributes, use the RCopy functions. If you know you only want to copy integer attributes, use the ICopy functions.

If the optional argument Vector is present and true, the vector architecture-friendly portions of this routine will be invoked.

If the optional argument sharedIndices is present, it should be the result of the call SharedIndices_(aVin, aVout, sharedIndices). Providing this argument speeds up this routine substantially. For example, you can compute a sharedIndices structure once for a given pair of AttrVects, then use that same structure for all copies between those two AttrVects (although note that a different sharedIndices variable would be needed if aVin and aVout were reversed). Note, however, that sharedIndices is ignored if either rList or iList are given.

N.B.: This routine will fail if the aVout is not initialized or if any of the specified attributes are not present in either aVout or aVin.


INTERFACE:

 
  subroutine Copy_(aVin, aVout, rList, TrList, iList, TiList, vector, sharedIndices)
USES:
       use m_die ,          only : die, warn
       use m_stdio ,        only : stderr
 
       use m_List,          only : GetSharedListIndices
       use m_List,          only : GetIndices => get_indices
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),             intent(in)    :: aVin
       character(len=*), optional, intent(in)    :: iList
       character(len=*), optional, intent(in)    :: rList
       character(len=*), optional, intent(in)    :: TiList
       character(len=*), optional, intent(in)    :: TrList
       logical, optional,          intent(in)    :: vector 
       type(AVSharedIndices), optional, intent(in) :: sharedIndices
OUTPUT PARAMETERS:
 
       type(AttrVect),             intent(inout) :: aVout
REVISION HISTORY:
   12Jun02 - R. Jacob <[email protected]> - initial version.
   13Jun02 - R. Jacob <[email protected]> - copy shared attributes
             if no attribute lists are specified.
   30Sep02 - R. Jacob <[email protected]> - new argument order with all
             optional arguments last
   19Feb02 - E. Ong <[email protected]> - new implementation using 
             new list function get_indices and faster memory copy  
   28Oct03 - R. Jacob <[email protected]> - add optional vector
             argument to use vector-friendly code provided by Fujitsu
   16Aug06 - R. Jacob <[email protected]> - split into 4 routines:
             RCopy_,RCopyL_,ICopy_,ICopyL_
   28Apr11 - W.J. Sacks <[email protected]> - added sharedIndices argument

2.1.25 Sort_ - Use Attributes as Keys to Generate an Index Permutation

The subroutine Sort_() uses a list of keys defined by the List key_list, searches for the appropriate integer or real attributes referenced by the items in key_list ( that is, it identifies the appropriate entries in aV%iList and aV%rList), and then uses these keys to generate a permutation perm that will put the entries of the attribute vector aV in lexicographic order as defined by key_list (the ordering in key_list being from left to right.

N.B.: This routine will fail if aV%iList and aV%rList share one or more common entries.

N.B.: This routine will fail if one of the sorting keys presented is not present in aV%iList nor aV%rList.


INTERFACE:

 
  subroutine Sort_(aV, key_list, perm, descend, perrWith, dieWith)
USES:
       use m_String,        only : String
       use m_String,        only : String_tochar => tochar
       use m_String,        only : String_clean => clean
       use m_List ,         only : List_allocated => allocated
       use m_List ,         only : List_index => index
       use m_List ,         only : List_nitem => nitem
       use m_List ,         only : List_get   => get
       use m_die ,          only : die
       use m_stdio ,        only : stderr
       use m_SortingTools , only : IndexSet
       use m_SortingTools , only : IndexSort
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),                  intent(in) :: aV
       type(List),                      intent(in) :: key_list
       logical, dimension(:), optional, intent(in) :: descend
       character(len=*),      optional, intent(in) :: perrWith
       character(len=*),      optional, intent(in) :: dieWith
OUTPUT PARAMETERS:
       integer, dimension(:),           pointer    :: perm
REVISION HISTORY:
   20Oct00 - J.W. Larson <[email protected]> - initial prototype
   25Apr01 - R.L. Jacob <[email protected]> - add -1 to make a
             backwards loop go backwards
   14Jun01 - J. Larson / E. Ong -- Fixed logic bug in REAL attribute
             sort (discovered by E. Ong), and cleaned up error / 
             shutdown logic.

2.1.26 Permute_ - Permute AttrVect Elements

The subroutine Permute_() uses a a permutation perm (which can be generated by the routine Sort_() in this module) to rearrange the entries in the attribute integer and real storage areas of the input attribute vector aV-aV%iAttr and aV%rAttr, respectively.


INTERFACE:

 
  subroutine Permute_(aV, perm, perrWith, dieWith)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
       use m_SortingTools , only : Permute
 
       implicit none
INPUT PARAMETERS:
       integer, dimension(:),           intent(in)    :: perm
       character(len=*),      optional, intent(in)    :: perrWith
       character(len=*),      optional, intent(in)    :: dieWith
INPUT/OUTPUT PARAMETERS:
       type(AttrVect),                  intent(inout) :: aV
REVISION HISTORY:
   23Oct00 - J.W. Larson <[email protected]> - initial prototype

2.1.27 Unpermute_ - Unpermute AttrVect Elements

The subroutine Unpermute_() uses a a permutation perm (which can be generated by the routine Sort_() in this module) to rearrange the entries in the attribute integer and real storage areas of the input attribute vector aV-aV%iAttr and aV%rAttr, respectively. This is meant to be called on an aV that has already been permuted but it could also be used to perform the inverse operation implied by perm on an unpermuted aV.


INTERFACE:

 
  subroutine Unpermute_(aV, perm, perrWith, dieWith)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
       use m_SortingTools , only : Unpermute
 
       implicit none
INPUT PARAMETERS:
       integer, dimension(:),           intent(in)    :: perm
       character(len=*),      optional, intent(in)    :: perrWith
       character(len=*),      optional, intent(in)    :: dieWith
INPUT/OUTPUT PARAMETERS:
       type(AttrVect),                  intent(inout) :: aV
REVISION HISTORY:
   23Nov05 - R. Jacob <[email protected]> - based on Permute

2.1.28 SortPermute_ - In-place Lexicographic Sort of an AttrVect

The subroutine SortPermute_() uses the routine Sort_() to create an index permutation perm that will place the AttrVect entries in the lexicographic order defined by the keys in the List variable key_list. This permutation is then used by the routine Permute_() to place the AttreVect entries in lexicographic order.


INTERFACE:

 
  subroutine SortPermute_(aV, key_list, descend, perrWith, dieWith)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       implicit none
INPUT PARAMETERS:
       type(List),                       intent(in)    :: key_list
       logical , dimension(:), optional, intent(in)    :: descend
       character(len=*),       optional, intent(in)    :: perrWith
       character(len=*),       optional, intent(in)    :: dieWith
INPUT/OUTPUT PARAMETERS:
       type(AttrVect),                   intent(inout) :: aV
REVISION HISTORY:
   24Oct00 - J.W. Larson <[email protected]> - initial prototype

2.1.29 aVaVSharedAttrIndexList_ - AttrVect shared attributes.

aVaVSharedAttrIndexList_() takes a pair of user-supplied AttrVect variables aV1 and aV2, and for choice of either REAL or INTEGER attributes (as specified literally in the input CHARACTER argument attrib) returns the number of shared attributes NumShared, and arrays of indices Indices1 and Indices2 to their storage locations in aV1 and aV2, respectively.

N.B.: This routine returns two allocated arrays--Indices1(:) and Indices2(:)--which must be deallocated once the user no longer needs them. Failure to do this will create a memory leak.


INTERFACE:

 
  subroutine aVaVSharedAttrIndexList_(aV1, aV2, attrib, NumShared, &
                                      Indices1, Indices2)
USES:
       use m_stdio
       use m_die,      only : MP_perr_die, die, warn
 
       use m_List,     only : GetSharedListIndices
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),        intent(in)  :: aV1   
       type(AttrVect),        intent(in)  :: aV2
       character(len=*),      intent(in)  :: attrib
OUTPUT PARAMETERS:
       integer,               intent(out) :: NumShared
       integer, dimension(:), pointer     :: Indices1
       integer, dimension(:), pointer     :: Indices2
REVISION HISTORY:
   07Feb01 - J.W. Larson <[email protected]> - initial version

2.1.30 SharedIndices_ - AttrVect shared attributes and auxiliary information

SharedIndices_() takes a pair of user-supplied AttrVect variables aV1 and aV2, and returns a structure of type AVSharedIndices (sharedIndices). This structure contains arrays of indices to the locations of the shared attributes, as well as auxiliary information. The structure contains information on both the REAL and INTEGER attributes. See documentation for the SharedIndicesOneType_ subroutine for some additional details, as much of the work is done there.

N.B.: The returned structure, sharedIndices, contains allocated arrays that must be deallocated once the user no longer needs them. This should be done through a call to cleanSharedIndices_.


INTERFACE:

 
  subroutine SharedIndices_(aV1, aV2, sharedIndices)
 
     implicit none
INPUT PARAMETERS:
       type(AttrVect),        intent(in)  :: aV1   
       type(AttrVect),        intent(in)  :: aV2
INPUT/OUTPUT PARAMETERS:
       type(AVSharedIndices), intent(inout) :: sharedIndices
REVISION HISTORY:
   28Apr11 - W.J. Sacks <[email protected]> - initial version

2.1.31 SharedIndicesOneType_ - AttrVect shared attributes and auxiliary information, for one data type

SharedIndicesOneType_() takes a pair of user-supplied AttrVect variables aV1 and aV2, and for choice of either REAL or INTEGER attributes (as specified literally in the input CHARACTER argument attrib) returns a structure of type AVSharedIndicesOneType ( sharedIndices). This structure contains arrays of indices to the locations of the shared attributes of the given type, as well as auxiliary information.

The aVindices1 and aVindices2 components of sharedIndices will be indices into aV1 and aV2, respectively.

N.B.: The returned structure, sharedIndices, contains allocated arrays that must be deallocated once the user no longer needs them. This should be done through a call to cleanSharedIndicesOneType_. Even if there are no attributes in common between aV1 and aV2, sharedIndices will still be initialized, and memory will still be allocated. Furthermore, if an already-initialized sharedIndices variable is to be given new values, cleanSharedIndicesOneType_ must be called before SharedIndicesOneType_ is called a second time, in order to prevent a memory leak.


INTERFACE:

 
  subroutine SharedIndicesOneType_(aV1, aV2, attrib, sharedIndices)
 
     implicit none
INPUT PARAMETERS:
       type(AttrVect),        intent(in)  :: aV1   
       type(AttrVect),        intent(in)  :: aV2
       character(len=*),      intent(in)  :: attrib
INPUT/OUTPUT PARAMETERS:
       type(AVSharedIndicesOneType), intent(inout) :: sharedIndices
REVISION HISTORY:
   28Apr11 - W.J. Sacks <[email protected]> - initial version

2.1.32 cleanSharedIndices_ - Deallocate allocated memory structures of an AVSharedIndices structure

This routine deallocates the allocated memory structures of the input/output AVSharedIndicesOneType argument sharedIndices, if they are currently associated. It also resets other components of this structure to a default state. The success (failure) of this operation is signified by a zero (non-zero) value of the optional INTEGER output argument stat. If clean_() is invoked without supplying stat, and any of the deallocation operations fail, the routine will terminate with an error message. If multiple errors occur, stat will give the error condition for the last error.


INTERFACE:

 
  subroutine cleanSharedIndices_(sharedIndices, stat)
 
     implicit none
INPUT/OUTPUT PARAMETERS:
       type(AVSharedIndices), intent(inout) :: sharedIndices
OUTPUT PARAMETERS:
       integer, optional, intent(out)       :: stat
REVISION HISTORY:
   28Apr11 - W.J. Sacks <[email protected]> - initial version

2.1.33 cleanSharedIndicesOneType_ - Deallocate allocated memory structures of an AVSharedIndicesOneType structure

This routine deallocates the allocated memory structures of the input/output AVSharedIndices argument sharedIndices, if they are currently associated. It also resets other components of this structure to a default state. The success (failure) of this operation is signified by a zero (non-zero) value of the optional INTEGER output argument stat. If clean_() is invoked without supplying stat, and any of the deallocation operations fail, the routine will terminate with an error message. If multiple errors occur, stat will give the error condition for the last error.


INTERFACE:

 
  subroutine cleanSharedIndicesOneType_(sharedIndices, stat)
USES:
       use m_die,      only : die
 
     implicit none
INPUT/OUTPUT PARAMETERS:
       type(AVSharedIndicesOneType), intent(inout) :: sharedIndices
OUTPUT PARAMETERS:
       integer, optional, intent(out)       :: stat
REVISION HISTORY:
   28Apr11 - W.J. Sacks <[email protected]> - initial version


2.2 Module m_AttrVectComms - MPI Communications Methods for the AttrVect (Source File: m_AttrVectComms.F90)

This module defines the communications methods for the AttrVect datatype (see the module m_AttrVect for more information about this class and its methods). MCT's communications are implemented in terms of the Message Passing Interface (MPI) standard, and we have as best as possible, made the interfaces to these routines appear as similar as possible to the corresponding MPI routines. For the AttrVect, we supply blocking point-to-point send and receive operations. We also supply the following collective operations: broadcast, gather, and scatter. The gather and scatter operations rely on domain decomposition descriptors that are defined elsewhere in MCT: the GlobalMap, which is a one-dimensional decomposition (see the MCT module m_GlobalMap for more details); and the GlobalSegMap, which is a segmented decomposition capable of supporting multidimensional domain decompositions (see the MCT module m_GlobalSegMap for more details).


INTERFACE:

  module m_AttrVectComms
USES:
       use m_AttrVect ! AttrVect class and its methods
 
       implicit none
 
       private	! except
 
       public :: gather		! gather all local vectors to the root
       public :: scatter		! scatter from the root to all PEs
       public :: bcast		! bcast from root to all PEs
       public :: send		! send an AttrVect
       public :: recv		! receive an AttrVect
 
     interface gather ; module procedure &
 	      GM_gather_, &
 	      GSM_gather_ 
     end interface
     interface scatter ; module procedure &
 	      GM_scatter_, &
 	      GSM_scatter_ 
     end interface
     interface bcast  ; module procedure bcast_  ; end interface
     interface send  ; module procedure send_  ; end interface
     interface recv  ; module procedure recv_  ; end interface
REVISION HISTORY:
   27Oct00 - J.W. Larson <[email protected]> - relocated routines
             from m_AttrVect to create this module.
   15Jan01 - J.W. Larson <[email protected]> - Added APIs for 
             GSM_gather_() and GSM_scatter_().
    9May01 - J.W. Larson <[email protected]> - Modified GM_scatter_
             so its communication model agrees with MPI_scatter().
             Also tidied up prologues in all module routines.
    7Jun01 - J.W. Larson <[email protected]> - Added send() 
             and recv().
    3Aug01 - E.T. Ong <[email protected]> - in GSM_scatter, call 
             GlobalMap_init with actual shaped array to satisfy
             Fortran 90 standard. See comment in subroutine.
   23Aug01 - E.T. Ong <[email protected]> - replaced assignment(=)
             with copy for list type to avoid compiler bugs in pgf90.
             Added more error checking in gsm scatter. Fixed minor bugs 
            in gsm and gm gather.
   13Dec01 - E.T. Ong <[email protected]> - GSM_scatter, allow users
             to scatter with a haloed GSMap. Fixed some bugs in 
             GM_scatter.
   19Dec01 - E.T. Ong <[email protected]> - allow bcast of an AttrVect
             with only an integer or real attribute.
   27Mar02 - J.W. Larson <[email protected]> - Corrected usage of
             m_die routines throughout this module.

2.2.1 send_ - Point-to-point Send of an AttrVect

This routine takes an input AttrVect argument inAV and sends it to processor dest on the communicator associated with the Fortran INTEGER MPI communicator handle comm. The overalll message is tagged by the input INTEGER argument TagBase. The success (failure) of this operation is reported in the zero (nonzero) optional output argument status.

N.B.: One must avoid assigning elsewhere the MPI tag values between TagBase and TagBase+7, inclusive. This is because send_() performs the send of the AttrVect as a series of eight send operations.


INTERFACE:

 
  subroutine send_(inAV, dest, TagBase, comm, status)
USES:
       use m_stdio
       use m_mpif90
       use m_die
 
       use m_List, only : List
       use m_List, only : List_allocated => allocated
       use m_List, only : List_nitem => nitem
       use m_List, only : List_send => send
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),     intent(in)  :: inAV
       integer,            intent(in)  :: dest
       integer,            intent(in)  :: TagBase
       integer,            intent(in)  :: comm
OUTPUT PARAMETERS:
       integer, optional,  intent(out) :: status
REVISION HISTORY:
    7Jun01 - J.W. Larson - initial version.
   13Jun01 - J.W. Larson <[email protected]> - Initialize status
             (if present).

2.2.2 recv_ - Point-to-point Receive of an AttrVect

This routine receives the output AttrVect argument outAV from processor source on the communicator associated with the Fortran INTEGER MPI communicator handle comm. The overall message is tagged by the input INTEGER argument TagBase. The success (failure) of this operation is reported in the zero (nonzero) optional output argument status.

N.B.: One must avoid assigning elsewhere the MPI tag values between TagBase and TagBase+7, inclusive. This is because recv_() performs the receive of the AttrVect as a series of eight receive operations.


INTERFACE:

 
  subroutine recv_(outAV, dest, TagBase, comm, status)
USES:
       use m_stdio
       use m_mpif90
       use m_die
 
       use m_List, only : List
       use m_List, only : List_nitem => nitem
       use m_List, only : List_recv => recv
 
       use m_AttrVect, only : AttrVect
 
       implicit none
INPUT PARAMETERS:
       integer,            intent(in)  :: dest
       integer,            intent(in)  :: TagBase
       integer,            intent(in)  :: comm
OUTPUT PARAMETERS:
       type(AttrVect),     intent(out) :: outAV
       integer, optional,  intent(out) :: status
REVISION HISTORY:
    7Jun01 - J.W. Larson - initial working version.
   13Jun01 - J.W. Larson <[email protected]> - Initialize status
             (if present).

2.2.3 GM_gather_ - Gather an AttrVect Distributed by a GlobalMap

This routine gathers a distributed AttrVect iV to the root process, and returns it in the output AttrVect argument oV. The decomposition of iV is described by the input GlobalMap argument GMap. The input INTEGER argument comm is the Fortran integer MPI communicator handle. The success (failure) of this operation corresponds to a zero (nonzero) value of the optional output INTEGER argument stat.


INTERFACE:

 
  subroutine GM_gather_(iV, oV, GMap, root, comm, stat)
USES:
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : FP
       use m_GlobalMap, only : GlobalMap
       use m_GlobalMap, only : GlobalMap_lsize => lsize
       use m_GlobalMap, only : GlobalMap_gsize => gsize
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero 
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nIAttr => nIAttr
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_clean => clean
       use m_FcComms,  only : fc_gatherv_int, fc_gatherv_fp
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),           intent(in)  :: iV
       type(GlobalMap),          intent(in)  :: GMap
       integer,                  intent(in)  :: root
       integer,                  intent(in)  :: comm
OUTPUT PARAMETERS:
       type(AttrVect),           intent(out) :: oV
       integer,        optional, intent(out) :: stat
REVISION HISTORY:
   15Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   27Oct00 - J.W. Larson <[email protected]> - relocated from
             m_AttrVect
   15Jan01 - J.W. Larson <[email protected]> - renamed GM_gather_
    9May01 - J.W. Larson <[email protected]> - tidied up prologue
   18May01 - R.L. Jacob <[email protected]> - use MP_Type function
             to determine type for mpi_gatherv
   31Jan09 - P.H. Worley <[email protected]> - replaced call to
             MPI_gatherv with call to flow controlled gather routines

2.2.4 GSM_gather_ - Gather an AttrVect Distributed by a GlobalSegMap

The routine GSM_gather_() takes a distributed input AttrVect argument iV, whose decomposition is described by the input GlobalSegMap argument GSMap, and gathers it to the output AttrVect argument oV. The gathered AttrVect oV is valid only on the root process specified by the input argument root. The communicator used to gather the data is specified by the argument comm. The success (failure) is reported in the zero (non-zero) value of the output argument stat.

GSM_gather_() converts the problem of gathering data according to a GlobalSegMap into the simpler problem of gathering data as specified by a GlobalMap. The GlobalMap variable GMap is created based on the local storage requirements for each distributed piece of iV. On the root, a complete (including halo points) gathered copy of iV is collected into the temporary AttrVect variable workV (the length of workV is the larger of GlobalSegMap_GlobalStorage(GSMap) or GlobalSegMap_GlobalSize(GSMap)). The variable workV is segmented by process, and segments are copied into it by process, but ordered in the same order the segments appear in GSMap. Once workV is loaded, the data are copied segment-by-segment to their appropriate locations in the output AttrVect oV.


INTERFACE:

 
  subroutine GSM_gather_(iV, oV, GSMap, root, comm, stat, rdefault, idefault)
USES:
   Message-passing environment utilities (mpeu) modules:
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only: FP
   GlobalSegMap and associated services:
       use m_GlobalSegMap, only : GlobalSegMap
       use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_id
       use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
       use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
       use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
       use m_GlobalSegMap, only : GlobalSegMap_haloed => haloed
       use m_GlobalSegMap, only : GlobalSegMap_GlobalStorage => GlobalStorage
   AttrVect and associated services:
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nIAttr => nIAttr
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_clean => clean 
   GlobalMap and associated services:
       use m_GlobalMap, only : GlobalMap
       use m_GlobalMap, only : GlobalMap_init => init
       use m_GlobalMap, only : GlobalMap_clean => clean
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),            intent(in)  :: iV
       type(GlobalSegMap),        intent(in)  :: GSMap
       integer,                   intent(in)  :: root
       integer,                   intent(in)  :: comm
       real(FP), optional,        intent(in)  :: rdefault
       integer, optional,         intent(in)  :: idefault
OUTPUT PARAMETERS:
       type(AttrVect),            intent(out) :: oV
       integer,        optional,  intent(out) :: stat
REVISION HISTORY:
   15Jan01 - J.W. Larson <[email protected]> - API specification.
   25Feb01 - J.W. Larson <[email protected]> - Prototype code.
   26Apr01 - R.L. Jacob <[email protected]> - add use statement for
             AttVect_clean
    9May01 - J.W. Larson <[email protected]> - tidied up prologue
   13Jun01 - J.W. Larson <[email protected]> - Initialize stat
             (if present).
   20Aug01 - E.T. Ong <[email protected]> - Added error checking for
             matching processors in gsmap and comm. Corrected
             current_pos assignment.
   23Nov01 - R. Jacob <[email protected]> - zero the oV before copying in
             gathered data.
   27Jul07 - R. Loy <[email protected]> - add Tony's suggested improvement
             for a default value in the output AV
   11Aug08 - R. Jacob <[email protected]> - add Pat Worley's faster way
             to initialize lns

2.2.5 GM_scatter_ - Scatter an AttrVect Using a GlobalMap

The routine GM_scatter_ takes an input AttrVect type iV (valid only on the root), and scatters it to a distributed AttrVect oV. The input GlobalMap argument GMap dictates how iV is scattered to oV. The success (failure) of this routine is reported in the zero (non-zero) value of the output argument stat.

N.B.: The output AttrVect argument oV represents dynamically allocated memory. When it is no longer needed, it should be deallocated by invoking AttrVect_clean() (see the module m_AttrVect for more details).


INTERFACE:

 
  subroutine GM_scatter_(iV, oV, GMap, root, comm, stat)
USES:
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : FP
 
       use m_List, only : List
       use m_List, only : List_copy => copy
       use m_List, only : List_bcast => bcast
       use m_List, only : List_clean => clean
       use m_List, only : List_nullify => nullify
       use m_List, only : List_nitem => nitem
 
       use m_GlobalMap, only : GlobalMap
       use m_GlobalMap, only : GlobalMap_lsize => lsize
       use m_GlobalMap, only : GlobalMap_gsize => gsize
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nIAttr => nIAttr
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_clean => clean
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),           intent(in)  :: iV
       type(GlobalMap),          intent(in)  :: GMap
       integer,                  intent(in)  :: root
       integer,                  intent(in)  :: comm
OUTPUT PARAMETERS:
       type(AttrVect),           intent(out) :: oV
       integer,        optional, intent(out) :: stat
REVISION HISTORY:
   21Apr98 - Jing Guo <guo@thunder> - initial prototype/prolog/code
   27Oct00 - J.W. Larson <[email protected]> - relocated from
             m_AttrVect
   15Jan01 - J.W. Larson <[email protected]> - renamed  GM_scatter_
    8Feb01 - J.W. Larson <[email protected]> - add logic to prevent
             empty calls (i.e. no data in buffer) to MPI_SCATTERV()
   27Apr01 - R.L. Jacob <[email protected]> - small bug fix to
             integer attribute scatter
    9May01 - J.W. Larson <[email protected]> - Re-vamped comms model
             to reflect MPI comms model for the scatter.  Tidied up
             the prologue, too.
   18May01 - R.L. Jacob <[email protected]> - use MP_Type function
             to determine type for mpi_scatterv
    8Aug01 - E.T. Ong <[email protected]> - replace list assignment(=)
             with list copy to avoid compiler errors in pgf90.
   13Dec01 - E.T. Ong <[email protected]> - allow scatter with an
             AttrVect containing only an iList or rList.

2.2.6 GSM_scatter_ - Scatter an AttrVect using a GlobalSegMap

The routine GSM_scatter_ takes an input AttrVect type iV (valid only on the root), and scatters it to a distributed AttrVect oV. The input GlobalSegMap argument GSMap dictates how iV is scattered to oV. The success (failure) of this routine is reported in the zero (non-zero) value of the output argument stat.

GSM_scatter_() converts the problem of scattering data according to a GlobalSegMap into the simpler problem of scattering data as specified by a GlobalMap. The GlobalMap variable GMap is created based on the local storage requirements for each distributed piece of iV. On the root, a complete (including halo points) copy of iV is stored in the temporary AttrVect variable workV (the length of workV is GlobalSegMap_GlobalStorage(GSMap)). The variable workV is segmented by process, and segments are copied into it by process, but ordered in the same order the segments appear in GSMap. Once workV is loaded, the data are scattered to the output AttrVect oV by a call to the routine GM_scatter_() defined in this module, with workV and GMap as the input arguments.

N.B.: This algorithm assumes that memory access times are much shorter than message-passing transmission times.

N.B.: The output AttrVect argument oV represents dynamically allocated memory. When it is no longer needed, it should be deallocated by invoking AttrVect_clean() (see the module m_AttrVect for more details).


INTERFACE:

 
  subroutine GSM_scatter_(iV, oV, GSMap, root, comm, stat)
USES:
   Environment utilities from mpeu:
 
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_List, only : List_nullify => nullify
 
   GlobalSegMap and associated services:
       use m_GlobalSegMap, only : GlobalSegMap
       use m_GlobalSegMap, only : GlobalSegMap_comp_id => comp_id
       use m_GlobalSegMap, only : GlobalSegMap_ngseg => ngseg
       use m_GlobalSegMap, only : GlobalSegMap_lsize => lsize
       use m_GlobalSegMap, only : GlobalSegMap_gsize => gsize
       use m_GlobalSegMap, only : GlobalSegMap_GlobalStorage => GlobalStorage
   AttrVect and associated services:
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero
       use m_AttrVect,  only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nIAttr => nIAttr
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_clean => clean 
   GlobalMap and associated services:
       use m_GlobalMap, only : GlobalMap
       use m_GlobalMap, only : GlobalMap_init => init
       use m_GlobalMap, only : GlobalMap_clean => clean
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),            intent(in)  :: iV
       type(GlobalSegMap),        intent(in)  :: GSMap
       integer,                   intent(in)  :: root
       integer,                   intent(in)  :: comm
OUTPUT PARAMETERS:
       type(AttrVect),            intent(out) :: oV
       integer,        optional,  intent(out) :: stat
REVISION HISTORY:
   15Jan01 - J.W. Larson <[email protected]> - API specification.
    8Feb01 - J.W. Larson <[email protected]> - Initial code.
   25Feb01 - J.W. Larson <[email protected]> - Bug fix--replaced
             call to GlobalSegMap_lsize with call to the new fcn.
             GlobalSegMap_ProcessStorage().
   26Apr01 - R.L. Jacob <[email protected]> - add use statement for
             AttVect_clean
   26Apr01 - J.W. Larson <[email protected]> - bug fixes--data 
             misalignment in use of the GlobalMap to compute the
             memory map into workV, and initialization of workV
             on all processes.
    9May01 - J.W. Larson <[email protected]> - tidied up prologue
   15May01 - Larson / Jacob <[email protected]> - stopped initializing
             workV on off-root processes (no longer necessary).
   13Jun01 - J.W. Larson <[email protected]> - Initialize stat
             (if present).
   20Jun01 - J.W. Larson <[email protected]> - Fixed a subtle bug
             appearing on AIX regarding the fact workV is uninitial-
             ized on non-root processes.  This is fixed by nullifying
             all the pointers in workV for non-root processes.
   20Aug01 - E.T. Ong <[email protected]> - Added argument check
             for matching processors in gsmap and comm.
   13Dec01 - E.T. Ong <[email protected]> - got rid of restriction 
             GlobalStorage(GSMap)==AttrVect_lsize(AV) to allow for
             GSMap to be haloed.
   11Aug08 - R. Jacob <[email protected]> - remove call to ProcessStorage
             and replace with faster algorithm provided by Pat Worley

2.2.7 bcast_ - Broadcast an AttrVect

This routine takes an AttrVect argument aV (at input, valid on the root only), and broadcasts it to all the processes associated with the communicator handle comm. The success (failure) of this routine is reported in the zero (non-zero) value of the output argument stat.

N.B.: The output (on non-root processes) AttrVect argument aV represents dynamically allocated memory. When it is no longer needed, it should be deallocated by invoking AttrVect_clean() (see the module m_AttrVect for details).


INTERFACE:

 
  subroutine bcast_(aV, root, comm, stat)
USES:
       use m_stdio
       use m_die
       use m_mpif90
       use m_String, only : String,bcast,char,String_clean
       use m_String, only : String_bcast => bcast
       use m_List, only : List_get => get
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nIAttr => nIAttr
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
 
       implicit none
INPUT PARAMETERS:
       integer,                  intent(in)    :: root
       integer,                  intent(in)    :: comm
INPUT/OUTPUT PARAMETERS:
       type(AttrVect),           intent(inout) :: aV ! (IN) on the root, 
                                                     ! (OUT) elsewhere
OUTPUT PARAMETERS:
       integer,        optional, intent(out)   :: stat
REVISION HISTORY:
   27Apr98 - Jing Guo <guo@thunder> - initial prototype/prologue/code
   27Oct00 - J.W. Larson <[email protected]> - relocated from
             m_AttrVect
    9May01 - J.W. Larson <[email protected]> - tidied up prologue
   18May01 - R.L. Jacob <[email protected]> - use MP_Type function
             to determine type for bcast
   19Dec01 - E.T. Ong <[email protected]> - adjusted for case of AV with 
             only integer or real attribute


2.3 Module m_AttrVectReduce - Local/Distributed AttrVect Reduction Ops. (Source File: m_AttrVectReduce.F90)

This module provides routines to perform reductions on the AttrVect datatype. These reductions can either be the types of operations supported by MPI (currently, summation, minimum and maximum are available) that are applied either to all the attributes (both integer and real), or specific reductions applicable only to the real attributes of an AttrVect. This module provides services for both local (i.e., one address space) and global (distributed) reductions. The type of reduction is defined through use of one of the public data members of this module:

Value Action
AttrVectSUM Sum
AttrVectMIN Minimum
AttrVectMAX Maximum


INTERFACE:

 
  module m_AttrVectReduce
USES:
       No modules are used in the declaration section of this module.
 
       implicit none
 
       private	! except
PUBLIC MEMBER FUNCTIONS:
 
       public :: LocalReduce            ! Local reduction of all attributes
       public :: LocalReduceRAttr       ! Local reduction of REAL attributes
       public :: AllReduce              ! AllReduce for distributed AttrVect
       public :: GlobalReduce           ! Local Reduce followed by AllReduce
       public :: LocalWeightedSumRAttr  ! Local weighted sum of 
                                        ! REAL attributes
       public :: GlobalWeightedSumRAttr ! Global weighted sum of REAL 
                                        ! attributes for a distrubuted 
                                        ! AttrVect
 
     interface LocalReduce ; module procedure LocalReduce_ ; end interface
     interface LocalReduceRAttr
        module procedure LocalReduceRAttr_ 
     end interface
     interface AllReduce
        module procedure AllReduce_ 
     end interface
     interface GlobalReduce
        module procedure GlobalReduce_ 
     end interface
     interface LocalWeightedSumRAttr; module procedure &
        LocalWeightedSumRAttrSP_, &
        LocalWeightedSumRAttrDP_
     end interface
     interface GlobalWeightedSumRAttr; module procedure &
        GlobalWeightedSumRAttrSP_, &
        GlobalWeightedSumRAttrDP_
     end interface
PUBLIC DATA MEMBERS:
 
     public :: AttrVectSUM
     public :: AttrVectMIN
     public :: AttrVectMAX
 
     integer, parameter :: AttrVectSUM = 1
     integer, parameter :: AttrVectMIN = 2
     integer, parameter :: AttrVectMAX = 3
REVISION HISTORY:
    7May02 - J.W. Larson <[email protected]> - Created module 
             using routines originally prototyped in m_AttrVect.

2.3.1 LocalReduce_ - Local Reduction of INTEGER and REAL Attributes

The subroutine LocalReduce_() takes the input AttrVect argument inAV, and reduces each of its integer and real attributes, returning them in the output AttrVect argument outAV (which is created by this routine). The type of reduction is defined by the input INTEGER argument action. Allowed values for action are defined as public data members to this module, and are summarized below:


Value Action
AttrVectSUM Sum
AttrVectMIN Minimum
AttrVectMAX Maximum

N.B.: The output AttrVect argument outAV is allocated memory, and must be destroyed by invoking the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine LocalReduce_(inAV, outAV, action)
USES:
       use m_realkinds,     only : FP
       use m_die ,          only : die
       use m_stdio ,        only : stderr
       use m_AttrVect,      only : AttrVect
       use m_AttrVect,      only : AttrVect_init => init
       use m_AttrVect,      only : AttrVect_zero => zero
       use m_AttrVect,      only : AttrVect_nIAttr => nIAttr
       use m_AttrVect,      only : AttrVect_nRAttr => nRAttr
       use m_AttrVect,      only : AttrVect_lsize => lsize
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),  intent(IN)  :: inAV
       integer,         intent(IN)  :: action
OUTPUT PARAMETERS:
       type(AttrVect),  intent(OUT) :: outAV
REVISION HISTORY:
   16Apr02 - J.W. Larson <[email protected]> - initial prototype

2.3.2 LocalReduceRAttr_ - Local Reduction of REAL Attributes

The subroutine LocalReduceRAttr_() takes the input AttrVect argument inAV, and reduces each of its REAL attributes, returning them in the output AttrVect argument outAV (which is created by this routine). The type of reduction is defined by the input INTEGER argument action. Allowed values for action are defined as public data members to this module (see the declaration section of m_AttrVect, and are summarized below:


Value Action
AttrVectSUM Sum
AttrVectMIN Minimum
AttrVectMAX Maximum

N.B.: The output AttrVect argument outAV is allocated memory, and must be destroyed by invoking the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

  subroutine LocalReduceRAttr_(inAV, outAV, action)
USES:
       use m_realkinds,     only : FP
 
       use m_die ,          only : die
       use m_stdio ,        only : stderr
 
       use m_List,          only : List
       use m_List,          only : List_copy => copy
       use m_List,          only : List_exportToChar => exportToChar
       use m_List,          only : List_clean => clean
 
       use m_AttrVect,      only : AttrVect
       use m_AttrVect,      only : AttrVect_init => init
       use m_AttrVect,      only : AttrVect_zero => zero
       use m_AttrVect,      only : AttrVect_nIAttr => nIAttr
       use m_AttrVect,      only : AttrVect_nRAttr => nRAttr
       use m_AttrVect,      only : AttrVect_lsize => lsize
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),               intent(IN)  :: inAV
       integer,                      intent(IN)  :: action
OUTPUT PARAMETERS:
       type(AttrVect),               intent(OUT) :: outAV
REVISION HISTORY:
   16Apr02 - J.W. Larson <[email protected]> - initial prototype
    6May02 - J.W. Larson <[email protected]> - added optional
             argument weights(:)
    8May02 - J.W. Larson <[email protected]> - modified interface
             to return it to being a pure reduction operation.
    9May02 - J.W. Larson <[email protected]> - renamed from 
             LocalReduceReals_() to LocalReduceRAttr_() to make
             the name more consistent with other module procedure
             names in this module.

2.3.3 AllReduce_ - Reduction of INTEGER and REAL Attributes

The subroutine AllReduce_() takes the distributed input AttrVect argument inAV, and performs a global reduction of all its attributes across the MPI communicator associated with the Fortran90 INTEGER handle comm, and returns these reduced values to all processes in the AttrVect argument outAV (which is created by this routine). The reduction operation is specified by the user, and must have one of the values listed in the table below:

Value Action
AttrVectSUM Sum
AttrVectMIN Minimum
AttrVectMAX Maximum

N.B.: The output AttrVect argument outAV is allocated memory, and must be destroyed by invoking the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine AllReduce_(inAV, outAV, ReductionOp, comm, ierr)
USES:
       use m_die
       use m_stdio ,        only : stderr
       use m_mpif90
 
       use m_List,          only : List
       use m_List,          only : List_exportToChar => exportToChar
       use m_List,          only : List_allocated => allocated
 
       use m_AttrVect,      only : AttrVect
       use m_AttrVect,      only : AttrVect_init => init
       use m_AttrVect,      only : AttrVect_zero => zero
       use m_AttrVect,      only : AttrVect_lsize => lsize
       use m_AttrVect,      only : AttrVect_nIAttr => nIAttr
       use m_AttrVect,      only : AttrVect_nRAttr => nRAttr
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),               intent(IN)  :: inAV
       integer,                      intent(IN)  :: ReductionOp
       integer,                      intent(IN)  :: comm
OUTPUT PARAMETERS:
       type(AttrVect),               intent(OUT) :: outAV
       integer,        optional,     intent(OUT) :: ierr
REVISION HISTORY:
    8May02 - J.W. Larson <[email protected]> - initial version.
    9Jul02 - J.W. Larson <[email protected]> - slight modification;
             use List_allocated() to determine if there is attribute
             data to be reduced (this patch is to support the Sun
             F90 compiler).

2.3.4 GlobalReduce_ - Reduction of INTEGER and REAL Attributes

The subroutine GlobalReduce_() takes the distributed input AttrVect argument inAV, and performs a local reduction of all its integer and real attributes, followed by a an AllReduce of all the result of the local reduction across the MPI communicator associated with the Fortran90 INTEGER handle comm, and returns these reduced values to all processes in the AttrVect argument outAV (which is created by this routine). The reduction operation is specified by the user, and must have one of the values listed in the table below:

Value Action
AttrVectSUM Sum
AttrVectMIN Minimum
AttrVectMAX Maximum

N.B.: The output AttrVect argument outAV is allocated memory, and must be destroyed by invoking the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine GlobalReduce_(inAV, outAV, ReductionOp, comm, ierr)
USES:
       use m_die
       use m_stdio ,        only : stderr
       use m_mpif90
 
       use m_AttrVect,      only : AttrVect
       use m_AttrVect,      only : AttrVect_clean => clean
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),               intent(IN)  :: inAV
       integer,                      intent(IN)  :: ReductionOp
       integer,                      intent(IN)  :: comm
OUTPUT PARAMETERS:
       type(AttrVect),               intent(OUT) :: outAV
       integer,        optional,     intent(OUT) :: ierr
REVISION HISTORY:
    6May03 - J.W. Larson <[email protected]> - initial version.

2.3.5 LocalWeightedSumRAttrSP_ - Local Weighted Sum of REAL Attributes

The subroutine LocalWeightedSumRAttr_() takes the input AttrVect argument inAV, and performs a weighted sum of each of its REAL attributes, returning them in the output AttrVect argument outAV (which is created by this routine and will contain no integer attributes). The weights used for the summation are provided by the user in the input argument Weights(:). If the sum of the weights is desired, this can be returned as an attribute in outAV if the optional CHARACTER argument WeightSumAttr is provided (which will be concatenated onto the list of real attributes in inAV).

N.B.: The argument WeightSumAttr must not be identical to any of the real attribute names in inAV.

N.B.: The output AttrVect argument outAV is allocated memory, and must be destroyed by invoking the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

  subroutine LocalWeightedSumRAttrSP_(inAV, outAV, Weights, WeightSumAttr)
USES:
       use m_die ,          only : die
       use m_stdio ,        only : stderr
       use m_realkinds,     only : SP, FP
 
       use m_List,          only : List
       use m_List,          only : List_init => init
       use m_List,          only : List_clean => clean
       use m_List,          only : List_exportToChar => exportToChar
       use m_List,          only : List_concatenate => concatenate
 
       use m_AttrVect,      only : AttrVect
       use m_AttrVect,      only : AttrVect_init => init
       use m_AttrVect,      only : AttrVect_zero => zero
       use m_AttrVect,      only : AttrVect_nIAttr => nIAttr
       use m_AttrVect,      only : AttrVect_nRAttr => nRAttr
       use m_AttrVect,      only : AttrVect_lsize => lsize
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),               intent(IN)  :: inAV
       real(SP), dimension(:),       pointer     :: Weights
       character(len=*),   optional, intent(IN)  :: WeightSumAttr
OUTPUT PARAMETERS:
       type(AttrVect),               intent(OUT) :: outAV
REVISION HISTORY:
    8May02 - J.W. Larson <[email protected]> - initial version.
   14Jun02 - J.W. Larson <[email protected]> - bug fix regarding
             accumulation of weights when invoked with argument
             weightSumAttr.  Now works in MCT unit tester.

2.3.6 GlobalWeightedSumRAttrSP_ - Global Weighted Sum of REAL Attributes

The subroutine GlobalWeightedSumRAttr_() takes the distributed input AttrVect argument inAV, and performs a weighted global sum across the MPI communicator associated with the Fortran90 INTEGER handle comm of each of its REAL attributes, returning the sums to each process in the AttrVect argument outAV (which is created by this routine and will contain no integer attributes). The weights used for the summation are provided by the user in the input argument weights(:). If the sum of the weights is desired, this can be returned as an attribute in outAV if the optional CHARACTER argument WeightSumAttr is provided (which will be concatenated onto the list of real attributes in inAV to form the list of real attributes for outAV).

N.B.: The argument WeightSumAttr must not be identical to any of the real attribute names in inAV.

N.B.: The output AttrVect argument outAV is allocated memory, and must be destroyed by invoking the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

  subroutine GlobalWeightedSumRAttrSP_(inAV, outAV, Weights, comm, &
                                     WeightSumAttr)
USES:
       use m_die
       use m_stdio ,        only : stderr
       use m_mpif90
       use m_realkinds,     only : SP
 
       use m_List,          only : List
       use m_List,          only : List_exportToChar => exportToChar
 
       use m_AttrVect,      only : AttrVect
       use m_AttrVect,      only : AttrVect_clean => clean
       use m_AttrVect,      only : AttrVect_lsize => lsize
 
       implicit none
INPUT PARAMETERS:
       type(AttrVect),               intent(IN)  :: inAV
       real(SP), dimension(:),       pointer     :: Weights
       integer,                      intent(IN)  :: comm
       character(len=*),   optional, intent(IN)  :: WeightSumAttr
OUTPUT PARAMETERS:
       type(AttrVect),               intent(OUT) :: outAV
REVISION HISTORY:
    8May02 - J.W. Larson <[email protected]> - initial version.



next up previous contents
Next: 3 Global Segment Map Up: 1 Basic API's and Previous: 1 MCTWorld   Contents
[email protected]