MOM6
mom_domains::pass_vector_complete Interface Reference

Detailed Description

Complete a halo update on a pair of arrays representing the two components of a vector.

Definition at line 74 of file MOM_domains.F90.

Private functions

subroutine pass_vector_complete_3d (id_update, u_cmpt, v_cmpt, MOM_dom, direction, stagger, halo, clock)
 pass_vector_complete_3d completes a halo update for a pair of three-dimensional arrays representing the compontents of a three-dimensional horizontal vector. More...
 
subroutine pass_vector_complete_2d (id_update, u_cmpt, v_cmpt, MOM_dom, direction, stagger, halo, clock)
 pass_vector_complete_2d completes a halo update for a pair of two-dimensional arrays representing the compontents of a two-dimensional horizontal vector. More...
 

Functions and subroutines

◆ pass_vector_complete_2d()

subroutine mom_domains::pass_vector_complete::pass_vector_complete_2d ( integer, intent(in)  id_update,
real, dimension(:,:), intent(inout)  u_cmpt,
real, dimension(:,:), intent(inout)  v_cmpt,
type(mom_domain_type), intent(inout)  MOM_dom,
integer, intent(in), optional  direction,
integer, intent(in), optional  stagger,
integer, intent(in), optional  halo,
integer, intent(in), optional  clock 
)
private

pass_vector_complete_2d completes a halo update for a pair of two-dimensional arrays representing the compontents of a two-dimensional horizontal vector.

Parameters
[in]id_updateThe integer id of this update which has been returned from a previous call to pass_var_start.
[in,out]u_cmptThe nominal zonal (u) component of the vector pair which is having its halos points exchanged.
[in,out]v_cmptThe nominal meridional (v) component of the vector pair which is having its halos points exchanged.
[in,out]mom_domThe MOM_domain_type containing the mpp_domain needed to determine where data should be sent.
[in]directionAn optional integer indicating which directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST, TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional scalars discretized at the typical vector component locations. For example, TO_EAST sends the data to the processor to the east, so the halos on the western side are filled. TO_ALL is the default if omitted.
[in]staggerAn optional flag, which may be one of A_GRID, BGRID_NE, or CGRID_NE, indicating where the two components of the vector are discretized. Omitting stagger is the same as setting it to CGRID_NE.
[in]haloThe size of the halo to update - the full halo by default.
[in]clockThe handle for a cpu time clock that should be started then stopped to time this routine.

Definition at line 781 of file MOM_domains.F90.

781  integer, intent(in) :: id_update !< The integer id of this update which has been
782  !! returned from a previous call to
783  !! pass_var_start.
784  real, dimension(:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
785  !! pair which is having its halos points
786  !! exchanged.
787  real, dimension(:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
788  !! vector pair which is having its halos points
789  !! exchanged.
790  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
791  !! needed to determine where data should be
792  !! sent.
793  integer, optional, intent(in) :: direction !< An optional integer indicating which
794  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
795  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
796  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
797  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
798  !! is the default if omitted.
799  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
800  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
801  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
802  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
803  !! halo by default.
804  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
805  !! started then stopped to time this routine.
806  ! Local variables
807  integer :: stagger_local
808  integer :: dirflag
809 
810  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
811 
812  stagger_local = cgrid_ne ! Default value for type of grid
813  if (present(stagger)) stagger_local = stagger
814 
815  dirflag = to_all ! 60
816  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
817 
818  if (present(halo) .and. mom_dom%thin_halo_updates) then
819  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
820  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local, &
821  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
822  else
823  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
824  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local)
825  endif
826 
827  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
828 

◆ pass_vector_complete_3d()

subroutine mom_domains::pass_vector_complete::pass_vector_complete_3d ( integer, intent(in)  id_update,
real, dimension(:,:,:), intent(inout)  u_cmpt,
real, dimension(:,:,:), intent(inout)  v_cmpt,
type(mom_domain_type), intent(inout)  MOM_dom,
integer, intent(in), optional  direction,
integer, intent(in), optional  stagger,
integer, intent(in), optional  halo,
integer, intent(in), optional  clock 
)
private

pass_vector_complete_3d completes a halo update for a pair of three-dimensional arrays representing the compontents of a three-dimensional horizontal vector.

Parameters
[in]id_updateThe integer id of this update which has been returned from a previous call to pass_var_start.
[in,out]u_cmptThe nominal zonal (u) component of the vector pair which is having its halos points exchanged.
[in,out]v_cmptThe nominal meridional (v) component of the vector pair which is having its halos points exchanged.
[in,out]mom_domThe MOM_domain_type containing the mpp_domain needed to determine where data should be sent.
[in]directionAn optional integer indicating which directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST, TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional scalars discretized at the typical vector component locations. For example, TO_EAST sends the data to the processor to the east, so the halos on the western side are filled. TO_ALL is the default if omitted.
[in]staggerAn optional flag, which may be one of A_GRID, BGRID_NE, or CGRID_NE, indicating where the two components of the vector are discretized. Omitting stagger is the same as setting it to CGRID_NE.
[in]haloThe size of the halo to update - the full halo by default.
[in]clockThe handle for a cpu time clock that should be started then stopped to time this routine.

Definition at line 835 of file MOM_domains.F90.

835  integer, intent(in) :: id_update !< The integer id of this update which has been
836  !! returned from a previous call to
837  !! pass_var_start.
838  real, dimension(:,:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
839  !! pair which is having its halos points
840  !! exchanged.
841  real, dimension(:,:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
842  !! vector pair which is having its halos points
843  !! exchanged.
844  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
845  !! needed to determine where data should be
846  !! sent.
847  integer, optional, intent(in) :: direction !< An optional integer indicating which
848  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
849  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
850  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
851  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
852  !! is the default if omitted.
853  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
854  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
855  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
856  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
857  !! halo by default.
858  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
859  !! started then stopped to time this routine.
860  ! Local variables
861  integer :: stagger_local
862  integer :: dirflag
863 
864  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
865 
866  stagger_local = cgrid_ne ! Default value for type of grid
867  if (present(stagger)) stagger_local = stagger
868 
869  dirflag = to_all ! 60
870  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
871 
872  if (present(halo) .and. mom_dom%thin_halo_updates) then
873  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
874  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local, &
875  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
876  else
877  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
878  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local)
879  endif
880 
881  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
882 

The documentation for this interface was generated from the following file: