Skip to content

Commit

Permalink
Add new FMS1 interfaces
Browse files Browse the repository at this point in the history
  • Loading branch information
MJHarrison-GFDL committed Jul 19, 2023
1 parent 0691b12 commit 3b1a7d7
Show file tree
Hide file tree
Showing 2 changed files with 63 additions and 2 deletions.
24 changes: 23 additions & 1 deletion config_src/infra/FMS1/MOM_coms_infra.F90
Original file line number Diff line number Diff line change
Expand Up @@ -25,7 +25,7 @@ module MOM_coms_infra
!> Communicate an array, string or scalar from one PE to others
interface broadcast
module procedure broadcast_char, broadcast_int32_0D, broadcast_int64_0D, broadcast_int1D
module procedure broadcast_real0D, broadcast_real1D, broadcast_real2D
module procedure broadcast_real0D, broadcast_real1D, broadcast_real2D, broadcast_real3D
end interface broadcast

!> Compute a checksum for a field distributed over a PE list. If no PE list is
Expand Down Expand Up @@ -260,6 +260,28 @@ subroutine broadcast_real2D(dat, length, from_PE, PElist, blocking)

end subroutine broadcast_real2D


!> Communicate a 3-D array of reals from one PE to others
subroutine broadcast_real3D(dat, length, from_PE, PElist, blocking)
real, dimension(:,:,:), intent(inout) :: dat !< The data to communicate and destination
integer, intent(in) :: length !< The total number of data elements
integer, optional, intent(in) :: from_PE !< The source PE, by default the root PE
integer, optional, intent(in) :: PElist(:) !< The list of participating PEs, by default the
!! active PE set as previously set via Set_PElist.
logical, optional, intent(in) :: blocking !< If true, barriers are added around the call

integer :: src_PE ! The processor that is sending the data
logical :: do_block ! If true add synchronizing barriers

do_block = .false. ; if (present(blocking)) do_block = blocking
if (present(from_PE)) then ; src_PE = from_PE ; else ; src_PE = root_PE() ; endif

if (do_block) call mpp_sync(PElist)
call mpp_broadcast(dat, length, src_PE, PElist)
if (do_block) call mpp_sync_self(PElist)

end subroutine broadcast_real3D

! field_chksum wrappers

!> Compute a checksum for a field distributed over a PE list. If no PE list is
Expand Down
41 changes: 40 additions & 1 deletion config_src/infra/FMS1/MOM_io_infra.F90
Original file line number Diff line number Diff line change
Expand Up @@ -57,7 +57,7 @@ module MOM_io_infra
!> Read a data field from a file
interface read_field
module procedure read_field_4d
module procedure read_field_3d
module procedure read_field_3d, read_field_3d_region
module procedure read_field_2d, read_field_2d_region
module procedure read_field_1d, read_field_1d_int
module procedure read_field_0d, read_field_0d_int
Expand Down Expand Up @@ -696,6 +696,45 @@ subroutine read_field_3d(filename, fieldname, data, MOM_Domain, &
endif ; endif
end subroutine read_field_3d

!> This routine uses the fms_io subroutine read_data to read a region from a distributed or
!! global 3-D data field named "fieldname" from file "filename".
subroutine read_field_3d_region(filename, fieldname, data, start, nread, MOM_domain, &
no_domain, scale)
character(len=*), intent(in) :: filename !< The name of the file to read
character(len=*), intent(in) :: fieldname !< The variable name of the data in the file
real, dimension(:,:,:), intent(inout) :: data !< The 3-dimensional array into which the data
!! should be read
integer, dimension(:), intent(in) :: start !< The starting index to read in each of 4
!! dimensions. For this 3-d read, the
!! 4th values are always 1.
integer, dimension(:), intent(in) :: nread !< The number of points to read in each of 4
!! dimensions. For this 3-d read, the
!! 4th values are always 1.
type(MOM_domain_type), &
optional, intent(in) :: MOM_Domain !< The MOM_Domain that describes the decomposition
logical, optional, intent(in) :: no_domain !< If present and true, this variable does not
!! use domain decomposion.
real, optional, intent(in) :: scale !< A scaling factor that the field is multiplied
!! by before it is returned.

if (present(MOM_Domain)) then
call read_data(filename, fieldname, data, start, nread, domain=MOM_Domain%mpp_domain, &
no_domain=no_domain)
else
call read_data(filename, fieldname, data, start, nread, no_domain=no_domain)
endif

if (present(scale)) then ; if (scale /= 1.0) then
if (present(MOM_Domain)) then
call rescale_comp_data(MOM_Domain, data, scale)
else
! Dangerously rescale the whole array
data(:,:,:) = scale*data(:,:,:)
endif
endif ; endif
end subroutine read_field_3d_region


!> This routine uses the fms_io subroutine read_data to read a distributed
!! 4-D data field named "fieldname" from file "filename". Valid values for
!! "position" include CORNER, CENTER, EAST_FACE and NORTH_FACE.
Expand Down

0 comments on commit 3b1a7d7

Please sign in to comment.