MOM6
mom_domains::pass_vector Interface Reference

Detailed Description

Definition at line 69 of file MOM_domains.F90.

Private functions

subroutine pass_vector_3d (u_cmpt, v_cmpt, MOM_dom, direction, stagger, complete, halo)
 
subroutine pass_vector_2d (u_cmpt, v_cmpt, MOM_dom, direction, stagger, complete, halo)
 

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

Definition at line 454 of file MOM_domains.F90.

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

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

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

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