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

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

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

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

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

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