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 79 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...
 

Detailed Description

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

Definition at line 79 of file MOM_domains.F90.

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 804 of file MOM_domains.F90.

804  integer, intent(in) :: id_update !< The integer id of this update which has been
805  !! returned from a previous call to
806  !! pass_var_start.
807  real, dimension(:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
808  !! pair which is having its halos points
809  !! exchanged.
810  real, dimension(:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
811  !! vector pair which is having its halos points
812  !! exchanged.
813  type(mom_domain_type), intent(inout) :: mom_dom !< The MOM_domain_type containing the mpp_domain
814  !! needed to determine where data should be
815  !! sent.
816  integer, optional, intent(in) :: direction !< An optional integer indicating which
817  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
818  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
819  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
820  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
821  !! is the default if omitted.
822  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
823  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
824  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
825  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
826  !! halo by default.
827  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
828  !! started then stopped to time this routine.
829  ! Local variables
830  integer :: stagger_local
831  integer :: dirflag
832 
833  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
834 
835  stagger_local = cgrid_ne ! Default value for type of grid
836  if (present(stagger)) stagger_local = stagger
837 
838  dirflag = to_all ! 60
839  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
840 
841  if (present(halo) .and. mom_dom%thin_halo_updates) then
842  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
843  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local, &
844  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
845  else
846  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
847  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local)
848  endif
849 
850  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
851 

◆ 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 858 of file MOM_domains.F90.

858  integer, intent(in) :: id_update !< The integer id of this update which has been
859  !! returned from a previous call to
860  !! pass_var_start.
861  real, dimension(:,:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
862  !! pair which is having its halos points
863  !! exchanged.
864  real, dimension(:,:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
865  !! vector pair which is having its halos points
866  !! exchanged.
867  type(mom_domain_type), intent(inout) :: mom_dom !< The MOM_domain_type containing the mpp_domain
868  !! needed to determine where data should be
869  !! sent.
870  integer, optional, intent(in) :: direction !< An optional integer indicating which
871  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
872  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
873  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
874  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
875  !! is the default if omitted.
876  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
877  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
878  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
879  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
880  !! halo by default.
881  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
882  !! started then stopped to time this routine.
883  ! Local variables
884  integer :: stagger_local
885  integer :: dirflag
886 
887  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
888 
889  stagger_local = cgrid_ne ! Default value for type of grid
890  if (present(stagger)) stagger_local = stagger
891 
892  dirflag = to_all ! 60
893  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
894 
895  if (present(halo) .and. mom_dom%thin_halo_updates) then
896  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
897  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local, &
898  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
899  else
900  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
901  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local)
902  endif
903 
904  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
905 

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