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

Detailed Description

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

Definition at line 59 of file MOM_domains.F90.

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

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

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

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

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