MOM6
mom_domains::pass_vector_complete Interface Reference

Detailed Description

Definition at line 85 of file MOM_domains.F90.

Private functions

subroutine pass_vector_complete_3d (id_update, u_cmpt, v_cmpt, MOM_dom, direction, stagger, halo)
 
subroutine pass_vector_complete_2d (id_update, u_cmpt, v_cmpt, MOM_dom, direction, stagger, halo)
 

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 
)
private
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.

Definition at line 821 of file MOM_domains.F90.

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

◆ 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 
)
private
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.

Definition at line 888 of file MOM_domains.F90.

888  integer, intent(in) :: id_update !< The integer id of this update which has been
889  !! returned from a previous call to
890  !! pass_var_start.
891  real, dimension(:,:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
892  !! pair which is having its halos points
893  !! exchanged.
894  real, dimension(:,:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
895  !! vector pair which is having its halos points
896  !! exchanged.
897  type(mom_domain_type), intent(inout) :: mom_dom !< The MOM_domain_type containing the mpp_domain
898  !! needed to determine where data should be
899  !! sent.
900  integer, optional, intent(in) :: direction !< An optional integer indicating which
901  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
902  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
903  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
904  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
905  !! is the default if omitted.
906  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
907  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
908  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
909  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
910  !! halo by default.
911 ! Arguments: id_update - The integer id of this update which has been returned
912 ! from a previous call to pass_var_start.
913 ! (inout) u_cmpt - The nominal zonal (u) component of the vector pair which
914 ! is having its halos points exchanged.
915 ! (inout) v_cmpt - The nominal meridional (v) component of the vector pair
916 ! which is having its halos points exchanged.
917 ! (in) MOM_dom - The MOM_domain_type containing the mpp_domain needed to
918 ! determine where data should be sent.
919 ! (in) direction - An optional integer indicating which directions the
920 ! data should be sent. It is TO_ALL or the sum of any of
921 ! TO_EAST, TO_WEST, TO_NORTH, and TO_SOUTH, possibly
922 ! plus SCALAR_PAIR if these are paired non-directional
923 ! scalars discretized at the typical vector component
924 ! locations. For example, TO_EAST sends the data to the
925 ! processor to the east, so the halos on the western
926 ! side are filled. TO_ALL is the default if omitted.
927 ! (in) stagger - An optional flag, which may be one of A_GRID, BGRID_NE,
928 ! or CGRID_NE, indicating where the two components of the
929 ! vector are discretized. Omitting stagger is the same as
930 ! setting it to CGRID_NE.
931 ! (in,opt) halo - The size of the halo to update - the full halo by default.
932 ! (return value) - The integer index for this update.
933  integer :: stagger_local
934  integer :: dirflag
935 
936  stagger_local = cgrid_ne ! Default value for type of grid
937  if (present(stagger)) stagger_local = stagger
938 
939  dirflag = to_all ! 60
940  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
941 
942  if (present(halo) .and. mom_dom%thin_halo_updates) then
943  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
944  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local, &
945  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
946  else
947  call mpp_complete_update_domains(id_update, u_cmpt, v_cmpt, &
948  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local)
949  endif
950 

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