MOM6
mom_domains::pass_vector_start Interface Reference

Detailed Description

Initiate 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

integer function pass_vector_start_3d (u_cmpt, v_cmpt, MOM_dom, direction, stagger, complete, halo, clock)
 pass_vector_start_3d starts a halo update for a pair of three-dimensional arrays representing the compontents of a three-dimensional horizontal vector. More...
 
integer function pass_vector_start_2d (u_cmpt, v_cmpt, MOM_dom, direction, stagger, complete, halo, clock)
 pass_vector_start_2d starts a halo update for a pair of two-dimensional arrays representing the compontents of a two-dimensional horizontal vector. More...
 

Detailed Description

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

Definition at line 74 of file MOM_domains.F90.

Functions and subroutines

◆ pass_vector_start_2d()

integer function mom_domains::pass_vector_start::pass_vector_start_2d ( 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,
logical, intent(in), optional  complete,
integer, intent(in), optional  halo,
integer, intent(in), optional  clock 
)
private

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

Parameters
[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]completeAn optional argument indicating whether the halo updates should be completed before progress resumes. Omitting complete is the same as setting complete to .true.
[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.
Returns
The integer index for this update.

Definition at line 691 of file MOM_domains.F90.

691  real, dimension(:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
692  !! pair which is having its halos points
693  !! exchanged.
694  real, dimension(:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
695  !! vector pair which is having its halos points
696  !! exchanged.
697  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
698  !! needed to determine where data should be
699  !! sent.
700  integer, optional, intent(in) :: direction !< An optional integer indicating which
701  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
702  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
703  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
704  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
705  !! is the default if omitted.
706  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
707  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
708  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
709  logical, optional, intent(in) :: complete !< An optional argument indicating whether the
710  !! halo updates should be completed before progress resumes.
711  !! Omitting complete is the same as setting complete to .true.
712  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
713  !! halo by default.
714  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
715  !! started then stopped to time this routine.
716  integer :: pass_vector_start_2d !< The integer index for this
717  !! update.
718 
719  ! Local variables
720  integer :: stagger_local
721  integer :: dirflag
722 
723  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
724 
725  stagger_local = cgrid_ne ! Default value for type of grid
726  if (present(stagger)) stagger_local = stagger
727 
728  dirflag = to_all ! 60
729  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
730 
731  if (present(halo) .and. mom_dom%thin_halo_updates) then
732  pass_vector_start_2d = mpp_start_update_domains(u_cmpt, v_cmpt, &
733  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local, &
734  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
735  else
736  pass_vector_start_2d = mpp_start_update_domains(u_cmpt, v_cmpt, &
737  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local)
738  endif
739 
740  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
741 

◆ pass_vector_start_3d()

integer function mom_domains::pass_vector_start::pass_vector_start_3d ( 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,
logical, intent(in), optional  complete,
integer, intent(in), optional  halo,
integer, intent(in), optional  clock 
)
private

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

Parameters
[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]completeAn optional argument indicating whether the halo updates should be completed before progress resumes. Omitting complete is the same as setting complete to .true.
[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.
Returns
The integer index for this update.

Definition at line 748 of file MOM_domains.F90.

748  real, dimension(:,:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
749  !! pair which is having its halos points
750  !! exchanged.
751  real, dimension(:,:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
752  !! vector pair which is having its halos points
753  !! exchanged.
754  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
755  !! needed to determine where data should be
756  !! sent.
757  integer, optional, intent(in) :: direction !< An optional integer indicating which
758  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
759  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
760  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
761  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
762  !! is the default if omitted.
763  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
764  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
765  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
766  logical, optional, intent(in) :: complete !< An optional argument indicating whether the
767  !! halo updates should be completed before progress resumes.
768  !! Omitting complete is the same as setting complete to .true.
769  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
770  !! halo by default.
771  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
772  !! started then stopped to time this routine.
773  integer :: pass_vector_start_3d !< The integer index for this
774  !! update.
775  ! Local variables
776  integer :: stagger_local
777  integer :: dirflag
778 
779  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
780 
781  stagger_local = cgrid_ne ! Default value for type of grid
782  if (present(stagger)) stagger_local = stagger
783 
784  dirflag = to_all ! 60
785  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
786 
787  if (present(halo) .and. mom_dom%thin_halo_updates) then
788  pass_vector_start_3d = mpp_start_update_domains(u_cmpt, v_cmpt, &
789  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local, &
790  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
791  else
792  pass_vector_start_3d = mpp_start_update_domains(u_cmpt, v_cmpt, &
793  mom_dom%mpp_domain, flags=dirflag, gridtype=stagger_local)
794  endif
795 
796  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
797 

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