next up previous contents
Next: 13 Merging of Flux Up: 2 High Level API's Previous: 11 Matrix Vector Multiplication   Contents

Subsections

12 Spatial Integration and Averaging

12.1 Module m_SpatialIntegral - Spatial Integrals and Averages using a GeneralGrid (Source File: m_SpatialIntegral.F90)

This module provides spatial integration and averaging services for the MCT. For a field $\Phi$ sampled at a point ${\bf x}$ in some multidimensional domain $\Omega$, the integral $I$ of $\Phi({\bf x})$ is

\begin{displaymath}I = \int_{\Omega} \Phi ({\bf x}) d\Omega .\end{displaymath}

The spatial average $A$ of $\Phi({\bf x})$ over $\Omega$ is

\begin{displaymath}A = {{ \int_{\Omega} \Phi ({\bf x}) d\Omega} \over
{ \int_{\Omega} d\Omega} }. \end{displaymath}

Since the AttrVect represents a discretized field, the integrals above are implemented as:

\begin{displaymath}I = \sum_{i=1}^N \Phi_i \Delta \Omega_i \end{displaymath}

and

\begin{displaymath}A = {{\sum_{i=1}^N \Phi_i \Delta \Omega_i } \over
{\sum_{i=1}^N \Delta \Omega_i } }, \end{displaymath}

where $N$ is the number of physical locations, $\Phi_i$ is the value of the field $\Phi$ at location $i$, and $\Delta \Omega_i$ is the spatial weight (lenghth element, cross-sectional area element, volume element, et cetera) at location $i$.

MCT extends the concept of integrals and area/volume averages to include masked integrals and averages. MCT recognizes both integer and real masks. An integer mask $M$ is a vector of integers (one corresponding to each physical location) with each element having value either zero or one. Integer masks are used to include/exclude data from averages or integrals. For example, if one were to compute globally averaged cloud amount over land (but not ocean nor sea-ice), one would assign a $1$ to each location on the land and a $0$ to each non-land location. A real mask $F$ is a vector of real numbers (one corresponding to each physical location) with each element having value within the closed interval $[0,1]$. .Real masks are used to represent fractional area/volume coverage at a location by a given component model. For example, if one wishes to compute area averages over sea-ice, one must include the ice fraction present at each point. Masked Integrals and averages are represented in the MCT by:

\begin{displaymath}I = \sum_{i=1}^N {\prod_{j=1}^J M_i} {\prod_{k=1}^K F_i}
\Phi_i \Delta \Omega_i \end{displaymath}

and

\begin{displaymath}A = {{\sum_{i=1}^N \bigg({\prod_{j=1}^J M_i}\bigg) \bigg( {\p...
...bigg) \bigg( {\prod_{k=1}^K F_i} \bigg)
\Delta \Omega_i } }, \end{displaymath}

where $J$ is the number of integer masks and $K$ is the number of real masks.

All of the routines in this module assume field data is stored in an attribute vector (AttrVect), and the integration/averaging is performed only on the REAL attributes. Physical coordinate grid and mask information is assumed to be stored as attributes in either a GeneralGrid, or pre-combined into a single integer mask and a single real mask.


INTERFACE:

 
  module m_SpatialIntegral
      
       implicit none
 
       private	! except
PUBLIC MEMBER FUNCTIONS:
 
       public :: SpatialIntegral        ! Spatial Integral
       public :: SpatialAverage         ! Spatial Area Average
 
       public :: MaskedSpatialIntegral  ! Masked Spatial Integral
       public :: MaskedSpatialAverage   ! MaskedSpatial Area Average
 
       public :: PairedSpatialIntegrals ! A Pair of Spatial 
                                       ! Integrals 
 
       public :: PairedSpatialAverages  ! A Pair of Spatial 
                                       ! Area Averages
 
       public :: PairedMaskedSpatialIntegrals ! A Pair of Masked
                                              ! Spatial Integrals 
 
       public :: PairedMaskedSpatialAverages ! A Pair of Masked
                                             ! Spatial Area Averages
 
       interface SpatialIntegral ; module procedure &
 	   SpatialIntegralRAttrGG_
       end interface
       interface SpatialAverage ; module procedure &
 	   SpatialAverageRAttrGG_ 
       end interface
       interface MaskedSpatialIntegral ; module procedure &
 	   MaskedSpatialIntegralRAttrGG_
       end interface
       interface MaskedSpatialAverage ; module procedure &
 	   MaskedSpatialAverageRAttrGG_
       end interface
       interface PairedSpatialIntegrals ; module procedure &
 	    PairedSpatialIntegralRAttrGG_
       end interface
       interface PairedSpatialAverages ; module procedure &
 	    PairedSpatialAverageRAttrGG_
       end interface
       interface PairedMaskedSpatialIntegrals ; module procedure &
 	    PairedMaskedIntegralRAttrGG_
       end interface
       interface PairedMaskedSpatialAverages ; module procedure &
 	    PairedMaskedAverageRAttrGG_
       end interface
REVISION HISTORY:
   	25Oct01 - J.W. Larson <[email protected]> - Initial version
          9May02 - J.W. Larson <[email protected]> - Massive Refactoring.
      10-14Jun02 - J.W. Larson <[email protected]> - Added Masked methods.
      17-18Jun02 - J.W. Larson <[email protected]> - Added Paired/Masked 
                   methods.
         18Jun02 - J.W. Larson <[email protected]> - Renamed module from 
                   m_GlobalIntegral to m_SpatialIntegral.
         15Jan03 - E.T. Ong <[email protected]> - Initialized real-only 
                   AttrVects using nullfied integer lists. This circuitous 
                   hack was required because the compaq compiler does not
                   compile the function AttrVectExportListToChar.

12.1.1 SpatialIntegralRAttrGG_ - Compute spatial integral.

This routine computes spatial integrals of the REAL attributes of the REAL attributes of the input AttrVect argument inAv. SpatialIntegralRAttrGG_() takes the input AttrVect argument inAv and computes the spatial integral using weights stored in the GeneralGrid argument GGrid and identified by the CHARACTER tag WeightTag. The integral of each REAL attribute is returned in the output AttrVect argument outAv. If SpatialIntegralRAttrGG_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the weights are also summed and stored in outAv (and can be referenced with the attribute tag defined by the argumentWeightTag. If SpatialIntegralRAttrGG_() is invoked with the optional INTEGER argument comm (a Fortran MPI communicator handle), the summation operations for the integral are completed on the local process, then reduced across the communicator, with all processes receiving the result.

N.B.: The local lengths of the AttrVect argument inAv and the GeneralGrid GGrid must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in GGrid.

N.B.: If SpatialIntegralRAttrGG_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the value of WeightTag must not conflict with any of the REAL attribute tags in inAv.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine SpatialIntegralRAttrGG_(inAv, outAv, GGrid, WeightTag, &
                                    SumWeights, comm)
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_realkinds, only : FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
 
       use m_GeneralGrid, only : GeneralGrid
       use m_GeneralGrid, only : GeneralGrid_lsize => lsize
       use m_GeneralGrid, only : GeneralGrid_indexRA => indexRA
       use m_GeneralGrid, only : GeneralGrid_exportRAttr => exportRAttr
 
       use m_SpatialIntegralV, only: SpatialIntegralV
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),    intent(IN) :: inAv
       type(GeneralGrid), intent(IN) :: GGrid
       character(len=*),  intent(IN) :: WeightTag
       logical, optional, intent(IN) :: SumWeights
       integer, optional, intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv
REVISION HISTORY:
   	06Feb02 - J.W. Larson <[email protected]> - initial version
   	09May02 - J.W. Larson <[email protected]> - Refactored and
                   renamed SpatialIntegralRAttrGG_().
         07Jun02 - J.W. Larson <[email protected]> - Bug fix and further
                   refactoring.

12.1.2 SpatialAverageRAttrGG_ - Compute spatial average.

This routine computes spatial averages of the REAL attributes of the input AttrVect argument inAv. SpatialAverageRAttrGG_() takes the input AttrVect argument inAv and computes the spatial average using weights stored in the GeneralGrid argument GGrid and identified by the CHARACTER tag WeightTag. The average of each REAL attribute is returned in the output AttrVect argument outAv. If SpatialAverageRAttrGG_() is invoked with the optional INTEGER argument comm (a Fortran MPI communicator handle), the summation operations for the average are completed on the local process, then reduced across the communicator, with all processes receiving the result.

N.B.: The local lengths of the AttrVect argument inAv and the GeneralGrid GGrid must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in GGrid.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine SpatialAverageRAttrGG_(inAv, outAv, GGrid, WeightTag, comm)
 
   ! USES:
 
       use m_realkinds, only : FP
    
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero 
       use m_AttrVect, only : AttrVect_clean => clean
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_indexRA => indexRA
 
       use m_GeneralGrid, only : GeneralGrid
 
       use m_List, only : List
       use m_List, only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),    intent(IN) :: inAv
       type(GeneralGrid), intent(IN) :: GGrid
       character(len=*),  intent(IN) :: WeightTag
       integer, optional, intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv
REVISION HISTORY:
   	08Feb02 - J.W. Larson <[email protected]> - initial version
   	08May02 - J.W. Larson <[email protected]> - minor modifications:
                   1) renamed the routine to GlobalAverageRAttrGG_
                   2) changed calls to reflect new routine name
                      GlobalIntegralRAttrGG_().
         18Jun02 - J.W. Larson <[email protected]> - Renamed routine to
                   SpatialAverageRAttrGG_().

12.1.3 MaskedSpatialIntegralRAttrGG_ - Masked spatial integral.

This routine computes masked spatial integrals of the REAL attributes of the input AttrVect argument inAv, returning the masked integrals in the output AttrVect outAv. All of the masking data are assumed stored in the input GeneralGrid argument GGrid. If integer masks are to be used, their integer attribute names in GGrid are named as a colon-delimited list in the optional CHARACTER input argument iMaskTags. Real masks (if desired) are referenced by their real attribute names in GGrid are named as a colon-delimited list in the optional CHARACTER input argument rMaskTags. The user specifies a choice of mask combination method with the input LOGICAL argument UseFastMethod. If ${\tt UseFastMethod} = {\tt .FALSE.}$ this routine checks each mask entry to ensure that the integer masks contain only ones and zeroes, and that entries in the real masks are all in the closed interval $[0,1]$. If ${\tt UseFastMethod} = {\tt .TRUE.}$, this routine performs direct products of the masks, assuming that the user has validated them in advance. The optional LOGICAL input argument SumWeights determines whether the masked sum of the spatial weights is computed and returned in outAv with the real attribute name supplied in the optional CHARACTER input argument WeightSumTag. This integral can either be a local (i.e. a global memory space operation), or a global distributed integral. The latter is the case if the optional input INTEGER argument comm is supplied (which corresponds to a Fortran MPI communicatior handle).

N.B.: The local lengths of the AttrVect argument inAv and the input GeneralGrid GGrid must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in GGrid.

N.B.: If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE.. In this case, the none of REAL attribute tags in inAv may be named the same as the string contained in WeightSumTag, which is an attribute name reserved for the sum of the weights in the output AttrVect outAv.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine MaskedSpatialIntegralRAttrGG_(inAv, outAv, GGrid, SpatialWeightTag, &
                                          iMaskTags, rMaskTags, UseFastMethod,  &
                                          SumWeights, WeightSumTag, comm)
 
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_realkinds, only : FP
 
       use m_String, only : String
       use m_String, only : String_toChar => toChar
       use m_String, only : String_clean => clean
 
       use m_List, only : List
       use m_List, only : List_init => init
       use m_List, only : List_clean => clean
       use m_List, only : List_nitem => nitem
       use m_List, only : List_get => get
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
 
       use m_GeneralGrid, only : GeneralGrid
       use m_GeneralGrid, only : GeneralGrid_lsize => lsize
       use m_GeneralGrid, only : GeneralGrid_indexRA => indexRA
       use m_GeneralGrid, only : GeneralGrid_exportIAttr => exportIAttr
       use m_GeneralGrid, only : GeneralGrid_exportRAttr => exportRAttr
 
       use m_AttrVectReduce, only : AttrVect_GlobalWeightedSumRAttr => &
 	                                         GlobalWeightedSumRAttr
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       use m_SpatialIntegralV, only : MaskedSpatialIntegralV
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),                  intent(IN) :: inAv
       type(GeneralGrid),               intent(IN) :: GGrid
       character(len=*),                intent(IN) :: SpatialWeightTag
       character(len=*),      optional, intent(IN) :: iMaskTags
       character(len=*),      optional, intent(IN) :: rMaskTags
       logical,                         intent(IN) :: UseFastMethod
       logical,               optional, intent(IN) :: SumWeights
       character(len=*),      optional, intent(IN) :: WeightSumTag
       integer,               optional, intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),                  intent(OUT) :: outAv
REVISION HISTORY:
   	11Jun02 - J.W. Larson <[email protected]> - initial version

12.1.4 MaskedSpatialAverageRAttrGG_ - Masked spatial average.

This routine computes masked spatial averages of the REAL attributes of the input AttrVect argument inAv, returning the masked averages in the output AttrVect outAv. All of the masking data are assumed stored in the input GeneralGrid argument GGrid. If integer masks are to be used, their integer attribute names in GGrid are named as a colon-delimited list in the optional CHARACTER input argument iMaskTags. Real masks (if desired) are referenced by their real attribute names in GGrid are named as a colon-delimited list in the optional CHARACTER input argument rMaskTags. The user specifies a choice of mask combination method with the input LOGICAL argument UseFastMethod. If ${\tt UseFastMethod} = {\tt .FALSE.}$ this routine checks each mask entry to ensure that the integer masks contain only ones and zeroes, and that entries in the real masks are all in the closed interval $[0,1]$. If ${\tt UseFastMethod} = {\tt .TRUE.}$, this routine performs direct products of the masks, assuming that the user has validated them in advance. This averaging can either be a local (equivalent to a global memory space operation), or a global distributed integral. The latter is the case if the optional input INTEGER argument comm is supplied (which corresponds to a Fortran MPI communicatior handle).

N.B.: The local lengths of the AttrVect argument inAv and the input GeneralGrid GGrid must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in GGrid.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine MaskedSpatialAverageRAttrGG_(inAv, outAv, GGrid, SpatialWeightTag, &
                                         iMaskTags, rMaskTags, UseFastMethod,  &
   				        comm)
 
   ! USES:
 
       use m_realkinds, only : FP
 
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero 
       use m_AttrVect, only : AttrVect_clean => clean
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_indexRA => indexRA
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
 
       use m_GeneralGrid, only : GeneralGrid
       use m_GeneralGrid, only : GeneralGrid_lsize => lsize
       use m_GeneralGrid, only : GeneralGrid_indexRA => indexRA
 
       use m_List, only : List
       use m_List, only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),                  intent(IN) :: inAv
       type(GeneralGrid),               intent(IN) :: GGrid
       character(len=*),                intent(IN) :: SpatialWeightTag
       character(len=*),      optional, intent(IN) :: iMaskTags
       character(len=*),      optional, intent(IN) :: rMaskTags
       logical,                         intent(IN) :: UseFastMethod
       integer,               optional, intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),                  intent(OUT) :: outAv
REVISION HISTORY:
   	12Jun02 - J.W. Larson <[email protected]> - initial version

12.1.5 PairedSpatialIntegralRAttrGG_ - Do two spatial integrals at once.

This routine computes spatial integrals of the REAL attributes of the REAL attributes of the input AttrVect arguments inAv1 and inAv2, returning the integrals in the output AttrVect arguments outAv1 and outAv2, respectively . The integrals of inAv1 and inAv2 are computed using spatial weights stored in the input GeneralGrid arguments GGrid1 and GGrid2, respectively. The spatial weights in in GGrid1 and GGrid2 are identified by the input CHARACTER arguments WeightTag1 and WeightTag2, respectively. If SpatialIntegralRAttrGG_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the weights are also summed and stored in outAv1 and outAv2, and can be referenced with the attribute tags defined by the arguments WeightTag1 and WeightTag2, respectively. This paired integral is implicitly a distributed operation (the whole motivation for pairing the integrals is to reduce communication latency costs), and the Fortran MPI communicator handle is defined by the input INTEGER argument comm. The summation is an AllReduce operation, with all processes receiving the global sum.

N.B.: The local lengths of the AttrVect argument inAv1 and the GeneralGrid GGrid1 must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv1 and the point weights stored in GGrid1. The same relationship must apply between inAv2 and GGrid2.

N.B.: If SpatialIntegralRAttrGG_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the value of WeightTag1 must not conflict with any of the REAL attribute tags in inAv1 and the value of WeightTag2 must not conflict with any of the REAL attribute tags in inAv2.

N.B.: The output AttrVect arguments outAv1 and outAv2 are allocated data structures. The user must deallocate them using the routine AttrVect_clean() when they are no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine PairedSpatialIntegralRAttrGG_(inAv1, outAv1, GGrid1, WeightTag1, &
                                          inAv2, outAv2, GGrid2, WeightTag2, &
 					 SumWeights, comm)
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_realkinds, only : FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
 
       use m_GeneralGrid, only : GeneralGrid
       use m_GeneralGrid, only : GeneralGrid_lsize => lsize
       use m_GeneralGrid, only : GeneralGrid_indexRA => indexRA
       use m_GeneralGrid, only : GeneralGrid_exportRAttr => exportRAttr
 
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       use m_SpatialIntegralV, only : PairedSpatialIntegralsV
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),    intent(IN) :: inAv1
       type(GeneralGrid), intent(IN) :: GGrid1
       character(len=*),  intent(IN) :: WeightTag1
       type(AttrVect),    intent(IN) :: inAv2
       type(GeneralGrid), intent(IN) :: GGrid2
       character(len=*),  intent(IN) :: WeightTag2
       logical, optional, intent(IN) :: SumWeights
       integer,           intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv1
       type(AttrVect),    intent(OUT) :: outAv2
REVISION HISTORY:
   	09May02 - J.W. Larson <[email protected]> - Initial version.
   	10Jun02 - J.W. Larson <[email protected]> - Refactored--now
                   built on top of PairedIntegralRAttrV_().

12.1.6 PairedSpatialAverageRAttrGG_ - Do two spatial averages at once.

This routine computes spatial averages of the REAL attributes of the REAL attributes of the input AttrVect arguments inAv1 and inAv2, returning the integrals in the output AttrVect arguments outAv1 and outAv2, respectively . The integrals of inAv1 and inAv2 are computed using spatial weights stored in the input GeneralGrid arguments GGrid1 and GGrid2, respectively. The spatial weights in in GGrid1 and GGrid2 are identified by the input CHARACTER arguments WeightTag1 and WeightTag2, respectively. This paired average is implicitly a distributed operation (the whole motivation for pairing the averages is to reduce communication latency costs), and the Fortran MPI communicator handle is defined by the input INTEGER argument comm. The summation is an AllReduce operation, with all processes receiving the global sum.

N.B.: The local lengths of the AttrVect argument inAv1 and the GeneralGrid GGrid1 must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv1 and the point weights stored in GGrid1. The same relationship must apply between inAv2 and GGrid2.

N.B.: The output AttrVect arguments outAv1 and outAv2 are allocated data structures. The user must deallocate them using the routine AttrVect_clean() when they are no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine PairedSpatialAverageRAttrGG_(inAv1, outAv1, GGrid1, WeightTag1, &
                                          inAv2, outAv2, GGrid2, WeightTag2, &
                                          comm)
   ! USES:
 
       use m_realkinds, only : FP
    
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero 
       use m_AttrVect, only : AttrVect_clean => clean
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_indexRA => indexRA
 
       use m_GeneralGrid, only : GeneralGrid
       use m_GeneralGrid, only : GeneralGrid_lsize => lsize
       use m_GeneralGrid, only : GeneralGrid_indexRA => indexRA
       use m_GeneralGrid, only : GeneralGrid_exportRAttr => exportRAttr
 
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       use m_List, only : List
       use m_List, only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),    intent(IN) :: inAv1
       type(GeneralGrid), intent(IN) :: GGrid1
       character(len=*),  intent(IN) :: WeightTag1
       type(AttrVect),    intent(IN) :: inAv2
       type(GeneralGrid), intent(IN) :: GGrid2
       character(len=*),  intent(IN) :: WeightTag2
       integer,           intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv1
       type(AttrVect),    intent(OUT) :: outAv2
REVISION HISTORY:
   	09May02 - J.W. Larson <[email protected]> - Initial version.
   	14Jun02 - J.W. Larson <[email protected]> - Bug fix to reflect
                   new interface to PairedSpatialIntegralRAttrGG_().

12.1.7 PairedMaskedIntegralRAttrGG_ - Do two masked integrals at once.

This routine computes a pair of masked spatial integrals of the REAL attributes of the input AttrVect arguments inAv and inAv2, returning the masked integrals in the output AttrVect outAv1 and outAv2, respectively. All of the spatial weighting and masking data for each set of integrals are assumed stored in the input GeneralGrid arguments GGrid and GGrid2. If integer masks are to be used, their integer attribute names in GGrid1 and GGrid2 are named as a colon-delimited lists in the optional CHARACTER input arguments iMaskTags1 and iMaskTags2, respectively. Real masks (if desired) are referenced by their real attribute names in GGrid1 and GGrid2 are named as colon-delimited lists in the optional CHARACTER input arguments rMaskTags1 and rMaskTags2, respectively. The user specifies a choice of mask combination method with the input LOGICAL argument UseFastMethod. If ${\tt UseFastMethod} = {\tt .FALSE.}$ this routine checks each mask entry to ensure that the integer masks contain only ones and zeroes, and that entries in the real masks are all in the closed interval $[0,1]$. If ${\tt UseFastMethod} = {\tt .TRUE.}$, this routine performs direct products of the masks, assuming that the user has validated them in advance. The optional LOGICAL input argument SumWeights determines whether the masked sum of the spatial weights is computed and returned in outAv1 and outAv2 with the real attribute names supplied in the CHARACTER input arguments SpatialWeightTag1, and SpatialWeightTag2, respectively. This paired integral is implicitly a distributed operation (the whole motivation for pairing the averages is to reduce communication latency costs), and the Fortran MPI communicator handle is defined by the input INTEGER argument comm. The summation is an AllReduce operation, with all processes receiving the global sum.

N.B.: The local lengths of the AttrVect argument inAv1 and the GeneralGrid GGrid1 must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv1 and the point weights stored in GGrid1. The same relationship must apply between inAv2 and GGrid2.

N.B.: If PairedMaskedIntegralRAttrGG_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the value of SpatialWeightTag1 must not conflict with any of the REAL attribute tags in inAv1 and the value of SpatialWeightTag2 must not conflict with any of the REAL attribute tags in inAv2.

N.B.: The output AttrVect arguments outAv1 and outAv2 are allocated data structures. The user must deallocate them using the routine AttrVect_clean() when they are no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine PairedMaskedIntegralRAttrGG_(inAv1, outAv1, GGrid1, &
                                          SpatialWeightTag1, rMaskTags1, &
                                          iMaskTags1, inAv2, outAv2, GGrid2, &
                                          SpatialWeightTag2, rMaskTags2, &
                                          iMaskTags2, UseFastMethod, &
                                          SumWeights, comm)
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_realkinds, only : FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
 
       use m_GeneralGrid, only : GeneralGrid
       use m_GeneralGrid, only : GeneralGrid_lsize => lsize
       use m_GeneralGrid, only : GeneralGrid_indexRA => indexRA
       use m_GeneralGrid, only : GeneralGrid_exportRAttr => exportRAttr
 
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),              intent(IN) :: inAv1
       type(GeneralGrid),           intent(IN) :: GGrid1
       character(len=*),            intent(IN) :: SpatialWeightTag1
       character(len=*),  optional, intent(IN) :: iMaskTags1
       character(len=*),  optional, intent(IN) :: rMaskTags1
       type(AttrVect),              intent(IN) :: inAv2
       type(GeneralGrid),           intent(IN) :: GGrid2
       character(len=*),            intent(IN) :: SpatialWeightTag2
       character(len=*),  optional, intent(IN) :: iMaskTags2
       character(len=*),  optional, intent(IN) :: rMaskTags2
       logical,                     intent(IN) :: UseFastMethod
       logical,           optional, intent(IN) :: SumWeights
       integer,                     intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv1
       type(AttrVect),    intent(OUT) :: outAv2
REVISION HISTORY:
   	17Jun02 - J.W. Larson <[email protected]> - Initial version.
   	19Jun02 - J.W. Larson <[email protected]> - Shortened the name
                   for compatibility with the Portland Group f90 compiler

12.1.8 PairedMaskedAverageRAttrGG_ - Do two masked averages at once.

This routine computes a pair of masked spatial averages of the REAL attributes of the input AttrVect arguments inAv and inAv2, returning the masked averagess in the output AttrVect outAv1 and outAv2, respectively. All of the spatial weighting and masking data for each set of averages are assumed stored in the input GeneralGrid arguments GGrid and GGrid2. If integer masks are to be used, their integer attribute names in GGrid1 and GGrid2 are named as a colon-delimited lists in the optional CHARACTER input arguments iMaskTags1 and iMaskTags2, respectively. Real masks (if desired) are referenced by their real attribute names in GGrid1 and GGrid2 are named as colon-delimited lists in the optional CHARACTER input arguments rMaskTags1 and rMaskTags2, respectively. The user specifies a choice of mask combination method with the input LOGICAL argument UseFastMethod. If ${\tt UseFastMethod} = {\tt .FALSE.}$ this routine checks each mask entry to ensure that the integer masks contain only ones and zeroes, and that entries in the real masks are all in the closed interval $[0,1]$. If ${\tt UseFastMethod} = {\tt .TRUE.}$, this routine performs direct products of the masks, assuming that the user has validated them in advance. This paired average is implicitly a distributed operation (the whole motivation for pairing the averages is to reduce communication latency costs), and the Fortran MPI communicator handle is defined by the input INTEGER argument comm. The summation is an AllReduce operation, with all processes receiving the global sum.

N.B.: The local lengths of the AttrVect argument inAv1 and the GeneralGrid GGrid1 must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv1 and the point weights stored in GGrid1. The same relationship must apply between inAv2 and GGrid2.

N.B.: The output AttrVect arguments outAv1 and outAv2 are allocated data structures. The user must deallocate them using the routine AttrVect_clean() when they are no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine PairedMaskedAverageRAttrGG_(inAv1, outAv1, GGrid1, &
                                         SpatialWeightTag1, rMaskTags1, &
                                         iMaskTags1, inAv2, outAv2, GGrid2, &
                                         SpatialWeightTag2, rMaskTags2, &
                                         iMaskTags2, UseFastMethod, &
                                         comm)
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
 
       use m_realkinds, only : FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero 
       use m_AttrVect, only : AttrVect_clean => clean
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
 
       use m_GeneralGrid, only : GeneralGrid
       use m_GeneralGrid, only : GeneralGrid_lsize => lsize
       use m_GeneralGrid, only : GeneralGrid_indexRA => indexRA
       use m_GeneralGrid, only : GeneralGrid_exportRAttr => exportRAttr
 
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       use m_List, only : List
       use m_List, only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),              intent(IN) :: inAv1
       type(GeneralGrid),           intent(IN) :: GGrid1
       character(len=*),            intent(IN) :: SpatialWeightTag1
       character(len=*),  optional, intent(IN) :: iMaskTags1
       character(len=*),  optional, intent(IN) :: rMaskTags1
       type(AttrVect),              intent(IN) :: inAv2
       type(GeneralGrid),           intent(IN) :: GGrid2
       character(len=*),            intent(IN) :: SpatialWeightTag2
       character(len=*),  optional, intent(IN) :: iMaskTags2
       character(len=*),  optional, intent(IN) :: rMaskTags2
       logical,                     intent(IN) :: UseFastMethod
       integer,                     intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv1
       type(AttrVect),    intent(OUT) :: outAv2
REVISION HISTORY:
   	17Jun02 - J.W. Larson <[email protected]> - Initial version.
   	19Jun02 - J.W. Larson <[email protected]> - Shortened the name
                   for compatibility with the Portland Group f90 compiler
   	25Jul02 - J.W. Larson E.T. Ong - Bug fix.  This routine was 
                   previously doing integrals rather than area averages.


12.2 Module m_SpatialIntegralV - Spatial Integrals and Averages using vectors of weights (Source File: m_SpatialIntegralV.F90)

This module provides spatial integration and averaging services for the MCT similar to those in m_SpatialIntegral except the weights are provided by an input vector instead of through a GeneralGrid. See the description for m_SpatialIntegral for more information

Paired masked spatial integrals and averages have not yet been implemented in vector form.


INTERFACE:

 
  module m_SpatialIntegralV
 
       implicit none
 
       private	! except
PUBLIC MEMBER FUNCTIONS:
 
       public :: SpatialIntegralV        ! Spatial Integral
       public :: SpatialAverageV         ! Spatial Area Average
 
       public :: MaskedSpatialIntegralV  ! Masked Spatial Integral
       public :: MaskedSpatialAverageV   ! MaskedSpatial Area Average
 
       public :: PairedSpatialIntegralsV ! A Pair of Spatial 
                                       ! Integrals 
 
       public :: PairedSpatialAveragesV  ! A Pair of Spatial 
                                       ! Area Averages
 
       interface SpatialIntegralV ; module procedure &
 	   SpatialIntegralRAttrVSP_, &
 	   SpatialIntegralRAttrVDP_
       end interface
       interface SpatialAverageV ; module procedure &
 	   SpatialAverageRAttrVSP_, &
 	   SpatialAverageRAttrVDP_
       end interface
       interface MaskedSpatialIntegralV ; module procedure &
 	   MaskedSpatialIntegralRAttrVSP_, &
 	   MaskedSpatialIntegralRAttrVDP_
       end interface
       interface MaskedSpatialAverageV ; module procedure &
 	   MaskedSpatialAverageRAttrVSP_, &
 	   MaskedSpatialAverageRAttrVDP_
       end interface
       interface PairedSpatialIntegralsV ; module procedure &
 	    PairedSpatialIntegralRAttrVSP_, &
 	    PairedSpatialIntegralRAttrVDP_
       end interface
       interface PairedSpatialAveragesV ; module procedure &
 	    PairedSpatialAverageRAttrVSP_, &
 	    PairedSpatialAverageRAttrVDP_
       end interface
REVISION HISTORY:
   	4Jan04 - R.Jacob <[email protected]> - move Vector versions of routines
                   from m_SpatialIntegral to this file.

12.2.1 SpatialIntegralRAttrVSP_ - Compute spatial integral.

This routine computes spatial integrals of the REAL attributes of the REAL attributes of the input AttrVect argument inAv. SpatialIntegralRAttrV_() takes the input AttrVect argument inAv and computes the spatial integral using weights stored in the input REAL array argument Weights. The integral of each REAL attribute is returned in the output AttrVect argument outAv. If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the weights are also summed and stored in outAv (and can be referenced with the attribute name WeightTag. If SpatialIntegralRAttrV_() is invoked with the optional INTEGER argument comm (a Fortran MPI communicator handle), the summation operations for the integral are completed on the local process, then reduced across the communicator, with all processes receiving the result.

N.B.: The local lengths of the AttrVect argument inAv and the input array Weights must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in Weights.

N.B.: If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE.. In this case, the none of REAL attribute tags in inAv may be named the same as the string contained in WeightTag, which is an attribute name reserved for the sum of the weights in the output AttrVect outAv.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine SpatialIntegralRAttrVSP_(inAv, outAv, Weights, SumWeights, &
                                   WeightTag, comm)
 
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : SP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
 
       use m_AttrVectReduce, only : AttrVect_GlobalWeightedSumRAttr => &
 	                                         GlobalWeightedSumRAttr
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),               intent(IN) :: inAv
       real(SP), dimension(:),       pointer    :: Weights
       logical,            optional, intent(IN) :: SumWeights
       character(len=*),   optional, intent(IN) :: WeightTag
       integer,            optional, intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv
REVISION HISTORY:
   	07Jun02 - J.W. Larson <[email protected]> - initial version

12.2.2 SpatialAverageRAttrVSP_ - Compute spatial average.

This routine computes spatial averages of the REAL attributes of the input AttrVect argument inAv. SpatialAverageRAttrV_() takes the input AttrVect argument inAv and computes the spatial average using weights stored in the REAL array Weights. The average of each REAL attribute is returned in the output AttrVect argument outAv. If SpatialAverageRAttrV_() is invoked with the optional INTEGER argument comm (a Fortran MPI communicator handle), the summation operations for the average are completed on the local process, then reduced across the communicator, with all processes receiving the result.

N.B.: The local lengths of the AttrVect argument inAv and the input array Weights must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in Weights.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine SpatialAverageRAttrVSP_(inAv, outAv, Weights, comm)
 
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : SP, FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero
       use m_AttrVect, only : AttrVect_clean => clean
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_indexRA => indexRA
 
       use m_List, only : List
       use m_List, only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),     intent(IN)  :: inAv
       real(SP), dimension(:), pointer :: Weights
       integer, optional,  intent(IN)  :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),     intent(OUT) :: outAv
REVISION HISTORY:
   	10Jun02 - J.W. Larson <[email protected]> - initial version

12.2.3 MaskedSpatialIntegralRAttrVSP_ - Masked spatial integral.

This routine computes masked spatial integrals of the REAL attributes of the input AttrVect argument inAv, returning the masked integrals in the output AttrVect argument outAv. The masked integral is computed using weights stored in the input REAL array argument SpatialWeights. Integer masking (if desired) is provided in the optional input INTEGER array iMask, and real masking (if desired) is provided in the optional input REAL array rMask. If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the weights are also summed and stored in outAv (and can be referenced with the attribute name defined by the optional input CHARACTER argument WeightSumTag. If SpatialIntegralRAttrV_() is invoked with the optional INTEGER argument comm (a Fortran MPI communicator handle), the summation operations for the integral are completed on the local process, then reduced across the communicator, with all processes receiving the result. Otherwise, the integral is assumed to be local (or equivalent to a global address space).

N.B.: The local lengths of the AttrVect argument inAv and the input array Weights must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in SpatialWeights.

N.B.: If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE.. In this case, the none of REAL attribute tags in inAv may be named the same as the string contained in WeightSumTag, which is an attribute name reserved for the sum of the weights in the output AttrVect outAv.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine MaskedSpatialIntegralRAttrVSP_(inAv, outAv, SpatialWeights, iMask, &
                                         rMask, UseFastMethod, SumWeights, &
                                         WeightSumTag, comm)
 
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : SP, FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
 
       use m_AttrVectReduce, only : AttrVect_GlobalWeightedSumRAttr => &
 	                                         GlobalWeightedSumRAttr
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),                  intent(IN) :: inAv
       real(SP),dimension(:),           pointer    :: SpatialWeights
       integer, dimension(:), optional, pointer    :: iMask
       real(SP),dimension(:), optional, pointer    :: rMask
       logical,                         intent(IN) :: UseFastMethod
       logical,               optional, intent(IN) :: SumWeights
       character(len=*),      optional, intent(IN) :: WeightSumTag
       integer,               optional, intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),                  intent(OUT) :: outAv
REVISION HISTORY:
   	10Jun02 - J.W. Larson <[email protected]> - initial version

12.2.4 MaskedSpatialAverageRAttrVSP_ - Masked spatial average.

[NEEDS **LOTS** of work...] This routine computes spatial integrals of the REAL attributes of the REAL attributes of the input AttrVect argument inAv. SpatialIntegralRAttrV_() takes the input AttrVect argument inAv and computes the spatial integral using weights stored in the input REAL array argument Weights. The integral of each REAL attribute is returned in the output AttrVect argument outAv. If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the weights are also summed and stored in outAv (and can be referenced with the attribute name WeightTag. If SpatialIntegralRAttrV_() is invoked with the optional INTEGER argument comm (a Fortran MPI communicator handle), the summation operations for the integral are completed on the local process, then reduced across the communicator, with all processes receiving the result.

N.B.: The local lengths of the AttrVect argument inAv and the input array Weights must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv and the point weights stored in Weights.

N.B.: If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE.. In this case, the none of REAL attribute tags in inAv may be named the same as the string contained in WeightTag, which is an attribute name reserved for the sum of the weights in the output AttrVect outAv.

N.B.: The output AttrVect argument outAv is an allocated data structure. The user must deallocate it using the routine AttrVect_clean() when it is no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine MaskedSpatialAverageRAttrVSP_(inAv, outAv, SpatialWeights, iMask, &
                                        rMask, UseFastMethod, comm)
 
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : SP, FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero
       use m_AttrVect, only : AttrVect_clean => clean
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_indexRA => indexRA
 
       use m_List, only : List
       use m_List, only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),                  intent(IN) :: inAv
       real(SP),  dimension(:),         pointer    :: SpatialWeights
       integer, dimension(:), optional, pointer    :: iMask
       real(SP),dimension(:), optional, pointer    :: rMask
       logical,                         intent(IN) :: UseFastMethod
       integer,               optional, intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),                  intent(OUT) :: outAv
REVISION HISTORY:
   	11Jun02 - J.W. Larson <[email protected]> - initial version

12.2.5 PairedSpatialIntegralRAttrVSP_ - Do two spatial integrals at once.

This routine computes spatial integrals of the REAL attributes of the REAL attributes of the input AttrVect arguments inAv1 and inAv2, returning the integrals in the output AttrVect arguments outAv1 and outAv2, respectively . The integrals of inAv1 and inAv2 are computed using spatial weights stored in the input REAL array arguments Weights1 and Weights2, respectively. If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the weights are also summed and stored in outAv1 and outAv2, and can be referenced with the attribute tags defined by the arguments WeightName1 and WeightName2, respectively. This paired integral is implicitly a distributed operation (the whole motivation for pairing the integrals is to reduce communication latency costs), and the Fortran MPI communicator handle is defined by the input INTEGER argument comm. The summation is an AllReduce operation, with all processes receiving the global sum.

N.B.: The local lengths of the AttrVect argument inAv1 and the input REAL array Weights1 must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv1 and the point weights stored in Weights. The same relationship must apply between inAv2 and Weights2.

N.B.: If SpatialIntegralRAttrV_() is invoked with the optional LOGICAL input argument SumWeights set as .TRUE., then the value of WeightName1 must not conflict with any of the REAL attribute tags in inAv1 and the value of WeightName2 must not conflict with any of the REAL attribute tags in inAv2.

N.B.: The output AttrVect arguments outAv1 and outAv2 are allocated data structures. The user must deallocate them using the routine AttrVect_clean() when they are no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine PairedSpatialIntegralRAttrVSP_(inAv1, outAv1, Weights1, WeightName1, &
                                         inAv2, outAv2, Weights2, WeightName2, &
                                         SumWeights, comm)
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : SP, FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
 
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),     intent(IN)  :: inAv1
       real(SP),dimension(:),pointer   :: Weights1
       character(len=*),   intent(IN)  :: WeightName1
       type(AttrVect),     intent(IN)  :: inAv2
       real(SP),dimension(:),pointer   :: Weights2
       character(len=*),   intent(IN)  :: WeightName2
       logical, optional,  intent(IN)  :: SumWeights
       integer,            intent(IN)  :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),     intent(OUT) :: outAv1
       type(AttrVect),     intent(OUT) :: outAv2
REVISION HISTORY:
   	10Jun02 - J.W. Larson <[email protected]> - Initial version.

12.2.6 PairedSpatialAverageRAttrVSP_ - Do two spatial averages at once.

This routine computes spatial averages of the REAL attributes of the REAL attributes of the input AttrVect arguments inAv1 and inAv2, returning the integrals in the output AttrVect arguments outAv1 and outAv2, respectively . The averages of inAv1 and inAv2 are computed using spatial weights stored in the input REAL array arguments Weights1 and Weights2, respectively. This paired average is implicitly a distributed operation (the whole motivation for pairing the integrals is to reduce communication latency costs), and the Fortran MPI communicator handle is defined by the input INTEGER argument comm. The summation is an AllReduce operation, with all processes receiving the global sum.

N.B.: The local lengths of the AttrVect argument inAv1 and the array Weights must be equal. That is, there must be a one-to-one correspondence between the field point values stored in inAv1 and the spatial weights stored in Weights

N.B.: The output AttrVect arguments outAv1 and outAv2 are allocated data structures. The user must deallocate them using the routine AttrVect_clean() when they are no longer needed. Failure to do so will result in a memory leak.


INTERFACE:

 
  subroutine PairedSpatialAverageRAttrVSP_(inAv1, outAv1, Weights1, inAv2, &
                                        outAv2, Weights2, comm)
   ! USES:
 
       use m_stdio
       use m_die
       use m_mpif90
       use m_realkinds, only : SP, FP
 
       use m_AttrVect, only : AttrVect
       use m_AttrVect, only : AttrVect_init => init
       use m_AttrVect, only : AttrVect_zero => zero
       use m_AttrVect, only : AttrVect_clean => clean
       use m_AttrVect, only : AttrVect_lsize => lsize
       use m_AttrVect, only : AttrVect_nRAttr => nRAttr
       use m_AttrVect, only : AttrVect_indexRA => indexRA
 
       use m_AttrVectReduce, only : AttrVect_LocalWeightedSumRAttr => &
 	                                         LocalWeightedSumRAttr
 
       use m_List, only : List
       use m_List, only : List_nullify => nullify
 
       implicit none
INPUT PARAMETERS:
 
       type(AttrVect),     intent(IN) :: inAv1
       real(SP),dimension(:),pointer  :: Weights1
       type(AttrVect),     intent(IN) :: inAv2
       real(SP),dimension(:),pointer  :: Weights2
       integer,            intent(IN) :: comm
OUTPUT PARAMETERS:
 
       type(AttrVect),    intent(OUT) :: outAv1
       type(AttrVect),    intent(OUT) :: outAv2
REVISION HISTORY:
   	09May02 - J.W. Larson <[email protected]> - Initial version.



next up previous contents
Next: 13 Merging of Flux Up: 2 High Level API's Previous: 11 Matrix Vector Multiplication   Contents
[email protected]