Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Use MPI_Bcast instead of multiple p2p messages to update nest from parent #272

Merged
merged 14 commits into from
Jan 22, 2024
Merged
Show file tree
Hide file tree
Changes from 12 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
20 changes: 5 additions & 15 deletions driver/fvGFS/atmosphere.F90
Original file line number Diff line number Diff line change
Expand Up @@ -165,7 +165,7 @@ module atmosphere_mod
mpp_npes, mpp_pe, mpp_chksum, &
mpp_get_current_pelist, &
mpp_set_current_pelist, &
mpp_sync, mpp_sync_self, mpp_send, mpp_recv
mpp_sync, mpp_sync_self, mpp_send, mpp_recv, mpp_broadcast
use mpp_parameter_mod, only: EUPDATE, WUPDATE, SUPDATE, NUPDATE
use mpp_domains_mod, only: CENTER, CORNER, NORTH, EAST, WEST, SOUTH
use mpp_domains_mod, only: domain2d, mpp_update_domains, mpp_global_field
Expand Down Expand Up @@ -2433,7 +2433,6 @@ subroutine fill_nested_grid_cpl(this_grid, proc_in)
logical, intent(in), optional :: proc_in

real, allocatable :: g_dat(:,:,:)
integer :: p, sending_proc
integer :: isd_p, ied_p, jsd_p, jed_p
integer :: isg, ieg, jsg, jeg
integer :: isc, iec, jsc, jec
Expand All @@ -2453,31 +2452,22 @@ subroutine fill_nested_grid_cpl(this_grid, proc_in)
allocate( g_dat(isg:ieg, jsg:jeg, 1) )

call timing_on('COMM_TOTAL')
sending_proc = Atm(this_grid)%parent_grid%pelist(1) + &
( Atm(this_grid)%neststruct%parent_tile-tile_fine(Atm(this_grid)%parent_grid%grid_number)+ &
Atm(this_grid)%parent_grid%flagstruct%ntiles-1 )*Atm(this_grid)%parent_grid%npes_per_tile
!if (Atm(this_grid)%neststruct%parent_proc .and. Atm(this_grid)%neststruct%parent_tile == Atm(this_grid)%parent_grid%global_tile) then
if (Atm(this_grid)%neststruct%parent_tile == Atm(this_grid)%parent_grid%global_tile) then
call mpp_global_field(Atm(this_grid)%parent_grid%domain, &
Atm(this_grid)%parent_grid%parent2nest_2d(isd_p:ied_p,jsd_p:jed_p), &
g_dat(isg:,jsg:,1), position=CENTER)
if (mpp_pe() == sending_proc) then
do p=1,size(Atm(this_grid)%pelist)
call mpp_send(g_dat, size(g_dat), Atm(this_grid)%pelist(p))
enddo
endif
endif
if (any(Atm(this_grid)%pelist == mpp_pe())) then
call mpp_recv(g_dat, size(g_dat), sending_proc)

if(Atm(this_grid)%BcastMember) then
call mpp_broadcast(g_dat, size(g_dat),Atm(this_grid)%Bcast_ranks(1), Atm(this_grid)%Bcast_ranks)
endif

call timing_off('COMM_TOTAL')
if (process) then
call fill_nested_grid(Atm(this_grid)%parent2nest_2d, g_dat(isg:,jsg:,1), &
Atm(this_grid)%neststruct%ind_h, Atm(this_grid)%neststruct%wt_h, &
0, 0, isg, ieg, jsg, jeg, Atm(this_grid)%bd)
endif

call mpp_sync_self
deallocate(g_dat)

end subroutine fill_nested_grid_cpl
Expand Down
8 changes: 7 additions & 1 deletion model/fv_arrays.F90
Original file line number Diff line number Diff line change
Expand Up @@ -1302,6 +1302,11 @@ module fv_arrays_mod

integer, allocatable, dimension(:) :: pelist

! These are set in fv_control_init() and used in fill_nested_grid_cpl()
! to replace numerous p2p MPI transfers with a single mpp_broadcast()
integer, allocatable :: Bcast_ranks(:)
logical :: BcastMember

type(fv_grid_bounds_type) :: bd

type(fv_regional_bc_bounds_type) :: regional_bc_bounds
Expand Down Expand Up @@ -1843,6 +1848,7 @@ end subroutine allocate_fv_atmos_type

!>@brief The subroutine 'deallocate_fv_atmos_type' deallocates the fv_atmos_type.
subroutine deallocate_fv_atmos_type(Atm)
use mpi
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Is this still needed?

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Not needed anymore since we are using mpp_broadcast() now instead of MPI_Bcast().

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Would like to see this use statement removed prior to accepting.

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Removed unneeded use of MPI module in model/fv_arrays.F90


implicit none
type(fv_atmos_type), intent(INOUT) :: Atm
Expand Down Expand Up @@ -2063,7 +2069,7 @@ subroutine deallocate_fv_atmos_type(Atm)
call deallocate_fv_nest_BC_type(Atm%neststruct%delz_BC)
endif
#endif

if(allocated(Atm%Bcast_ranks)) deallocate(Atm%Bcast_ranks)
end if

if (Atm%flagstruct%grid_type < 4) then
Expand Down
22 changes: 21 additions & 1 deletion model/fv_control.F90
Original file line number Diff line number Diff line change
Expand Up @@ -452,7 +452,6 @@ subroutine fv_control_init(Atm, dt_atmos, this_grid, grids_on_this_pe, p_split,
allocate(global_pelist(npes))
call mpp_get_current_pelist(global_pelist, commID=global_commID) ! for commID


allocate(grids_master_procs(ngrids))
pecounter = 0
allocate(grids_on_this_pe(ngrids))
Expand Down Expand Up @@ -734,6 +733,27 @@ subroutine fv_control_init(Atm, dt_atmos, this_grid, grids_on_this_pe, p_split,
Atm(this_grid)%neststruct%parent_proc = ANY(Atm(this_grid)%neststruct%child_grids) !ANY(tile_coarse == Atm(this_grid)%global_tile)
Atm(this_grid)%neststruct%child_proc = ASSOCIATED(Atm(this_grid)%parent_grid) !this means a nested grid

if(n>1) then
call mpp_set_current_pelist( global_pelist )
do n=2,ngrids
! Construct the MPI communicators that are used in fill_nested_grid_cpl()
allocate(Atm(n)%Bcast_ranks(1+size(Atm(n)%pelist)))
! parent grid sending rank within Bcast_ranks array
Atm(n)%Bcast_ranks(1)=Atm(n)%parent_grid%pelist(1) + &
( Atm(n)%neststruct%parent_tile-tile_fine(Atm(n)%parent_grid%grid_number)+ &
Atm(n)%parent_grid%flagstruct%ntiles-1 )*Atm(n)%parent_grid%npes_per_tile

bensonr marked this conversation as resolved.
Show resolved Hide resolved
Atm(n)%Bcast_ranks(2:(1+size(Atm(n)%pelist)))=Atm(n)%pelist ! Receivers
Atm(n)%BcastMember=.false.
if(any(mpp_pe() == Atm(n)%Bcast_ranks)) then
Atm(n)%BcastMember=.true.
endif

call mpp_declare_pelist(Atm(n)%Bcast_ranks(:))
enddo
call mpp_set_current_pelist(Atm(this_grid)%pelist)
endif

if (ngrids > 1) call setup_update_regions
if (Atm(this_grid)%neststruct%nestbctype > 1) then
call mpp_error(FATAL, 'nestbctype > 1 not yet implemented')
Expand Down
16 changes: 3 additions & 13 deletions model/molecular_diffusion.F90
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ module molecular_diffusion_mod
use constants_mod, only: rdgas, cp_air
use fv_mp_mod, only: is_master
use mpp_mod, only: stdlog, input_nml_file
use fms_mod, only: check_nml_error, open_namelist_file, close_file
use fms_mod, only: check_nml_error
use fv_grid_utils_mod, only: g_sum
use mpp_domains_mod, only: domain2d
use fv_arrays_mod, only: fv_grid_type, fv_grid_bounds_type
Expand Down Expand Up @@ -128,20 +128,10 @@ subroutine read_namelist_molecular_diffusion_nml(nml_filename,ncnst,nwat)

unit = stdlog()

#ifdef INTERNAL_FILE_NML

! Read molecular_diffusion namelist
read (input_nml_file,molecular_diffusion_nml,iostat=ios)
ierr = check_nml_error(ios,'molecular_diffusion_nml')
read (input_nml_file,molecular_diffusion_nml,iostat=ios)
ierr = check_nml_error(ios,'molecular_diffusion_nml')

#else
! Read molecular_diffusion namelist
f_unit = open_namelist_file(nml_filename)
rewind (f_unit)
read (f_unit,molecular_diffusion_nml,iostat=ios)
ierr = check_nml_error(ios,'molecular_diffusion_nml')
call close_file(f_unit)
#endif
write(unit, nml=molecular_diffusion_nml)
call molecular_diffusion_init(ncnst,nwat)

Expand Down
24 changes: 7 additions & 17 deletions model/multi_gases.F90
Original file line number Diff line number Diff line change
Expand Up @@ -42,7 +42,7 @@ module multi_gases_mod
use constants_mod, only: rdgas, rvgas, cp_air
use fv_mp_mod, only: is_master
use mpp_mod, only: stdlog, input_nml_file
use fms_mod, only: check_nml_error, open_namelist_file, close_file
use fms_mod, only: check_nml_error


implicit none
Expand Down Expand Up @@ -154,25 +154,15 @@ subroutine read_namelist_multi_gases_nml(nml_filename,ncnst,nwat)
ri(1) = rvgas
cpi(0) = cp_air
cpi(1) = 4*cp_air
#ifdef INTERNAL_FILE_NML

! Read multi_gases namelist
read (input_nml_file,multi_gases_nml,iostat=ios)
ierr = check_nml_error(ios,'multi_gases_nml')
! Read multi_gases namelist
read (input_nml_file,multi_gases_nml,iostat=ios)
ierr = check_nml_error(ios,'multi_gases_nml')

#else
! Read multi_gases namelist
f_unit = open_namelist_file(nml_filename)
write(unit, nml=multi_gases_nml)
call multi_gases_init(ncnst,nwat)

rewind (f_unit)
read (f_unit,multi_gases_nml,iostat=ios)
ierr = check_nml_error(ios,'multi_gases_nml')
call close_file(f_unit)
#endif
write(unit, nml=multi_gases_nml)
call multi_gases_init(ncnst,nwat)

return
return
end subroutine read_namelist_multi_gases_nml


Expand Down
8 changes: 6 additions & 2 deletions tools/fv_grid_tools.F90
Original file line number Diff line number Diff line change
Expand Up @@ -148,7 +148,7 @@ module fv_grid_tools_mod
use fms2_io_mod, only: file_exists, variable_exists, open_file, read_data, &
get_global_attribute, get_variable_attribute, &
close_file, get_mosaic_tile_grid, FmsNetcdfFile_t
use mosaic_mod, only: get_mosaic_ntiles
use mosaic2_mod, only: get_mosaic_ntiles

implicit none
private
Expand Down Expand Up @@ -227,14 +227,18 @@ subroutine read_grid(Atm, grid_file, ndims, nregions, ng)

call get_mosaic_tile_grid(atm_hgrid, atm_mosaic, Atm%domain)

if (open_file(Grid_input, atm_mosaic, "read")) then
ntiles = get_mosaic_ntiles(Grid_input)
call close_file(Grid_input)
endif

bensonr marked this conversation as resolved.
Show resolved Hide resolved
grid_form = "none"
if (open_file(Grid_input, atm_hgrid, "read")) then
call get_global_attribute(Grid_input, "history", attvalue)
if( index(attvalue, "gnomonic_ed") > 0) grid_form = "gnomonic_ed"
if(grid_form .NE. "gnomonic_ed") call mpp_error(FATAL, &
"fv_grid_tools(read_grid): the grid should be 'gnomonic_ed' when reading from grid file, contact developer")

ntiles = get_mosaic_ntiles(atm_mosaic)
if( .not. Atm%gridstruct%bounded_domain) then !<-- The regional setup has only 1 tile so do not shutdown in that case.
if(ntiles .NE. 6) call mpp_error(FATAL, &
'fv_grid_tools(read_grid): ntiles should be 6 in mosaic file '//trim(atm_mosaic) )
Expand Down
6 changes: 3 additions & 3 deletions tools/fv_iau_mod.F90
Original file line number Diff line number Diff line change
Expand Up @@ -37,7 +37,7 @@

module fv_iau_mod

use fms_mod, only: file_exist
use fms2_io_mod, only: file_exists
use mpp_mod, only: mpp_error, FATAL, NOTE, mpp_pe
use mpp_domains_mod, only: domain2d

Expand Down Expand Up @@ -186,7 +186,7 @@ subroutine IAU_initialize (IPD_Control, IAU_Data,Init_parm)
npz = IPD_Control%levs
fname = 'INPUT/'//trim(IPD_Control%iau_inc_files(1))

if( file_exist(fname) ) then
if( file_exists(fname) ) then
call open_ncfile( fname, ncid ) ! open the file
call get_ncdim1( ncid, 'lon', im)
call get_ncdim1( ncid, 'lat', jm)
Expand Down Expand Up @@ -457,7 +457,7 @@ subroutine read_iau_forcing(IPD_Control,increments,fname)

npz = IPD_Control%levs

if( file_exist(fname) ) then
if( file_exists(fname) ) then
call open_ncfile( fname, ncid ) ! open the file
else
call mpp_error(FATAL,'==> Error in read_iau_forcing: Expected file '&
Expand Down
2 changes: 0 additions & 2 deletions tools/fv_restart.F90
Original file line number Diff line number Diff line change
Expand Up @@ -168,7 +168,6 @@ module fv_restart_mod
use mpp_domains_mod, only: mpp_global_field
use fv_treat_da_inc_mod, only: read_da_inc
use fms2_io_mod, only: file_exists, set_filename_appendix, FmsNetcdfFile_t, open_file, close_file
use fms_io_mod, only: fmsset_filename_appendix=> set_filename_appendix
use coarse_grained_restart_files_mod, only: fv_io_write_restart_coarse
use fv_regional_mod, only: write_full_fields
#ifdef MULTI_GASES
Expand Down Expand Up @@ -289,7 +288,6 @@ subroutine fv_restart(fv_domain, Atm, seconds, days, cold_start, grid_type, this
if (Atm(n)%neststruct%nested .and. n==this_grid) then
write(gnn,'(A4, I2.2)') "nest", Atm(n)%grid_number
call set_filename_appendix(gnn)
call fmsset_filename_appendix(gnn)
endif

!3preN. Topography BCs for nest, including setup for blending
Expand Down
2 changes: 1 addition & 1 deletion tools/test_cases.F90
Original file line number Diff line number Diff line change
Expand Up @@ -117,7 +117,7 @@ module test_cases_mod

use mpp_mod, only: mpp_error, FATAL, mpp_root_pe, mpp_broadcast, mpp_sum
use mpp_mod, only: stdlog, input_nml_file
use fms_mod, only: check_nml_error, close_file, open_namelist_file
use fms_mod, only: check_nml_error
use mpp_domains_mod, only: mpp_update_domains, domain2d
use mpp_parameter_mod, only: AGRID_PARAM=>AGRID,CGRID_NE_PARAM=>CGRID_NE, &
SCALAR_PAIR
Expand Down