MOM6
mom_domains::pass_vector Interface Reference

Detailed Description

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

Definition at line 54 of file MOM_domains.F90.

Private functions

subroutine pass_vector_3d (u_cmpt, v_cmpt, MOM_dom, direction, stagger, complete, halo, clock)
 pass_vector_3d does a halo update for a pair of three-dimensional arrays representing the compontents of a three-dimensional horizontal vector. More...
 
subroutine pass_vector_2d (u_cmpt, v_cmpt, MOM_dom, direction, stagger, complete, halo, clock)
 pass_vector_2d does 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_2d()

subroutine mom_domains::pass_vector::pass_vector_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_2d does 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.

Definition at line 464 of file MOM_domains.F90.

464  real, dimension(:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
465  !! pair which is having its halos points
466  !! exchanged.
467  real, dimension(:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
468  !! vector pair which is having its halos points
469  !! exchanged.
470  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
471  !! needed to determine where data should be
472  !! sent.
473  integer, optional, intent(in) :: direction !< An optional integer indicating which
474  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
475  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
476  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
477  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
478  !! is the default if omitted.
479  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
480  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
481  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
482  logical, optional, intent(in) :: complete !< An optional argument indicating whether the
483  !! halo updates should be completed before progress resumes.
484  !! Omitting complete is the same as setting complete to .true.
485  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
486  !! halo by default.
487  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
488  !! started then stopped to time this routine.
489 
490  ! Local variables
491  integer :: stagger_local
492  integer :: dirflag
493  logical :: block_til_complete
494 
495  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
496 
497  stagger_local = cgrid_ne ! Default value for type of grid
498  if (present(stagger)) stagger_local = stagger
499 
500  dirflag = to_all ! 60
501  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
502  block_til_complete = .true.
503  if (present(complete)) block_til_complete = complete
504 
505  if (present(halo) .and. mom_dom%thin_halo_updates) then
506  call mpp_update_domains(u_cmpt, v_cmpt, mom_dom%mpp_domain, flags=dirflag, &
507  gridtype=stagger_local, complete = block_til_complete, &
508  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
509  else
510  call mpp_update_domains(u_cmpt, v_cmpt, mom_dom%mpp_domain, flags=dirflag, &
511  gridtype=stagger_local, complete = block_til_complete)
512  endif
513 
514  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
515 

◆ pass_vector_3d()

subroutine mom_domains::pass_vector::pass_vector_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_3d does 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.

Definition at line 610 of file MOM_domains.F90.

610  real, dimension(:,:,:), intent(inout) :: u_cmpt !< The nominal zonal (u) component of the vector
611  !! pair which is having its halos points
612  !! exchanged.
613  real, dimension(:,:,:), intent(inout) :: v_cmpt !< The nominal meridional (v) component of the
614  !! vector pair which is having its halos points
615  !! exchanged.
616  type(MOM_domain_type), intent(inout) :: MOM_dom !< The MOM_domain_type containing the mpp_domain
617  !! needed to determine where data should be
618  !! sent.
619  integer, optional, intent(in) :: direction !< An optional integer indicating which
620  !! directions the data should be sent. It is TO_ALL or the sum of any of TO_EAST, TO_WEST,
621  !! TO_NORTH, and TO_SOUTH, possibly plus SCALAR_PAIR if these are paired non-directional
622  !! scalars discretized at the typical vector component locations. For example, TO_EAST sends
623  !! the data to the processor to the east, so the halos on the western side are filled. TO_ALL
624  !! is the default if omitted.
625  integer, optional, intent(in) :: stagger !< An optional flag, which may be one of A_GRID,
626  !! BGRID_NE, or CGRID_NE, indicating where the two components of the vector are
627  !! discretized. Omitting stagger is the same as setting it to CGRID_NE.
628  logical, optional, intent(in) :: complete !< An optional argument indicating whether the
629  !! halo updates should be completed before progress resumes.
630  !! Omitting complete is the same as setting complete to .true.
631  integer, optional, intent(in) :: halo !< The size of the halo to update - the full
632  !! halo by default.
633  integer, optional, intent(in) :: clock !< The handle for a cpu time clock that should be
634  !! started then stopped to time this routine.
635 
636  ! Local variables
637  integer :: stagger_local
638  integer :: dirflag
639  logical :: block_til_complete
640 
641  if (present(clock)) then ; if (clock>0) call cpu_clock_begin(clock) ; endif
642 
643  stagger_local = cgrid_ne ! Default value for type of grid
644  if (present(stagger)) stagger_local = stagger
645 
646  dirflag = to_all ! 60
647  if (present(direction)) then ; if (direction > 0) dirflag = direction ; endif
648  block_til_complete = .true.
649  if (present(complete)) block_til_complete = complete
650 
651  if (present(halo) .and. mom_dom%thin_halo_updates) then
652  call mpp_update_domains(u_cmpt, v_cmpt, mom_dom%mpp_domain, flags=dirflag, &
653  gridtype=stagger_local, complete = block_til_complete, &
654  whalo=halo, ehalo=halo, shalo=halo, nhalo=halo)
655  else
656  call mpp_update_domains(u_cmpt, v_cmpt, mom_dom%mpp_domain, flags=dirflag, &
657  gridtype=stagger_local, complete = block_til_complete)
658  endif
659 
660  if (present(clock)) then ; if (clock>0) call cpu_clock_end(clock) ; endif
661 

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