Skip to content

Commit

Permalink
Allocate shelf_sfc_mass_flux for dynamic ice shelf configurations
Browse files Browse the repository at this point in the history
  • Loading branch information
MJHarrison-GFDL committed Jan 13, 2022
1 parent 6da5c9b commit c1a9131
Show file tree
Hide file tree
Showing 2 changed files with 34 additions and 4 deletions.
29 changes: 27 additions & 2 deletions src/core/MOM_forcing_type.F90
Original file line number Diff line number Diff line change
Expand Up @@ -170,6 +170,8 @@ module MOM_forcing_type
!! exactly 0 away from shelves or on land.
real, pointer, dimension(:,:) :: iceshelf_melt => NULL() !< Ice shelf melt rate (positive)
!! or freezing (negative) [R Z T-1 ~> kg m-2 s-1]
real, pointer, dimension(:,:) :: shelf_sfc_mass_flux => NULL() !< Ice shelf surface mass flux
!! deposition from the atmosphere. [R Z T-1 ~> kg m-2 s-1]

! Scalars set by surface forcing modules
real :: vPrecGlobalAdj = 0. !< adjustment to restoring vprec to zero out global net [kg m-2 s-1]
Expand Down Expand Up @@ -377,7 +379,7 @@ module MOM_forcing_type
! Iceberg + Ice shelf diagnostic handles
integer :: id_ustar_ice_cover = -1
integer :: id_frac_ice_cover = -1

integer :: id_shelf_sfc_mass_flux = -1
! wave forcing diagnostics handles.
integer :: id_lamult = -1
!>@}
Expand Down Expand Up @@ -2099,6 +2101,12 @@ subroutine fluxes_accumulate(flux_tmp, fluxes, G, wt2, forces)
fluxes%iceshelf_melt(i,j) = flux_tmp%iceshelf_melt(i,j)
enddo ; enddo
endif
if (associated(fluxes%shelf_sfc_mass_flux) &
.and. associated(flux_tmp%shelf_sfc_mass_flux)) then
do i=isd,ied ; do j=jsd,jed
fluxes%shelf_sfc_mass_flux(i,j) = flux_tmp%shelf_sfc_mass_flux(i,j)
enddo ; enddo
endif
if (associated(fluxes%frac_shelf_h) .and. associated(flux_tmp%frac_shelf_h)) then
do i=isd,ied ; do j=jsd,jed
fluxes%frac_shelf_h(i,j) = flux_tmp%frac_shelf_h(i,j)
Expand Down Expand Up @@ -2910,6 +2918,10 @@ subroutine forcing_diagnostics(fluxes_in, sfc_state, G_in, US, time_end, diag, h
if ((handles%id_ustar_ice_cover > 0) .and. associated(fluxes%ustar_shelf)) &
call post_data(handles%id_ustar_ice_cover, fluxes%ustar_shelf, diag)

if ((handles%id_shelf_sfc_mass_flux > 0) &
.and. associated(fluxes%shelf_sfc_mass_flux)) &
call post_data(handles%id_shelf_sfc_mass_flux, fluxes%shelf_sfc_mass_flux, diag)

! wave forcing ===============================================================
if (handles%id_lamult > 0) &
call post_data(handles%id_lamult, fluxes%lamult, diag)
Expand All @@ -2928,7 +2940,8 @@ end subroutine forcing_diagnostics

!> Conditionally allocate fields within the forcing type
subroutine allocate_forcing_by_group(G, fluxes, water, heat, ustar, press, &
shelf, iceberg, salt, fix_accum_bug, cfc, waves)
shelf, iceberg, salt, fix_accum_bug, cfc, waves, &
shelf_sfc_accumulation)
type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
type(forcing), intent(inout) :: fluxes !< A structure containing thermodynamic forcing fields
logical, optional, intent(in) :: water !< If present and true, allocate water fluxes
Expand All @@ -2942,14 +2955,21 @@ subroutine allocate_forcing_by_group(G, fluxes, water, heat, ustar, press, &
!! accumulation of ustar_gustless
logical, optional, intent(in) :: cfc !< If present and true, allocate cfc fluxes
logical, optional, intent(in) :: waves !< If present and true, allocate wave fields
logical, optional, intent(in) :: shelf_sfc_accumulation !< If present and true, and shelf is true,
!! then allocate surface flux deposition from the atmosphere
!! over ice shelves and ice sheets.

! Local variables
integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB
logical :: heat_water
logical :: shelf_sfc_acc

isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
IsdB = G%IsdB ; IedB = G%IedB ; JsdB = G%JsdB ; JedB = G%JedB

shelf_sfc_acc=.false.
if (present(shelf_sfc_accumulation)) shelf_sfc_acc=shelf_sfc_accumulation

call myAlloc(fluxes%ustar,isd,ied,jsd,jed, ustar)
call myAlloc(fluxes%ustar_gustless,isd,ied,jsd,jed, ustar)

Expand Down Expand Up @@ -2990,6 +3010,7 @@ subroutine allocate_forcing_by_group(G, fluxes, water, heat, ustar, press, &
call myAlloc(fluxes%frac_shelf_h,isd,ied,jsd,jed, shelf)
call myAlloc(fluxes%ustar_shelf,isd,ied,jsd,jed, shelf)
call myAlloc(fluxes%iceshelf_melt,isd,ied,jsd,jed, shelf)
call myAlloc(fluxes%shelf_sfc_mass_flux,isd,ied,jsd,jed, shelf_sfc_acc)

!These fields should only on allocated when iceberg area is being passed through the coupler.
call myAlloc(fluxes%ustar_berg,isd,ied,jsd,jed, iceberg)
Expand Down Expand Up @@ -3257,6 +3278,8 @@ subroutine deallocate_forcing_type(fluxes)
if (associated(fluxes%ustar_tidal)) deallocate(fluxes%ustar_tidal)
if (associated(fluxes%ustar_shelf)) deallocate(fluxes%ustar_shelf)
if (associated(fluxes%iceshelf_melt)) deallocate(fluxes%iceshelf_melt)
if (associated(fluxes%shelf_sfc_mass_flux)) &
deallocate(fluxes%shelf_sfc_mass_flux)
if (associated(fluxes%frac_shelf_h)) deallocate(fluxes%frac_shelf_h)
if (associated(fluxes%ustar_berg)) deallocate(fluxes%ustar_berg)
if (associated(fluxes%area_berg)) deallocate(fluxes%area_berg)
Expand Down Expand Up @@ -3355,6 +3378,7 @@ subroutine rotate_forcing(fluxes_in, fluxes, turns)
call rotate_array(fluxes_in%frac_shelf_h, turns, fluxes%frac_shelf_h)
call rotate_array(fluxes_in%ustar_shelf, turns, fluxes%ustar_shelf)
call rotate_array(fluxes_in%iceshelf_melt, turns, fluxes%iceshelf_melt)
call rotate_array(fluxes_in%shelf_sfc_mass_flux, turns, fluxes%shelf_sfc_mass_flux)
endif

if (do_iceberg) then
Expand Down Expand Up @@ -3600,6 +3624,7 @@ subroutine homogenize_forcing(fluxes, G)
call homogenize_field_t(fluxes%frac_shelf_h, G)
call homogenize_field_t(fluxes%ustar_shelf, G)
call homogenize_field_t(fluxes%iceshelf_melt, G)
call homogenize_field_t(fluxes%shelf_sfc_mass_flux, G)
endif

if (do_iceberg) then
Expand Down
9 changes: 7 additions & 2 deletions src/ice_shelf/MOM_ice_shelf.F90
Original file line number Diff line number Diff line change
Expand Up @@ -188,7 +188,8 @@ module MOM_ice_shelf
id_h_shelf = -1, id_h_mask = -1, &
id_surf_elev = -1, id_bathym = -1, &
id_area_shelf_h = -1, &
id_ustar_shelf = -1, id_shelf_mass = -1, id_mass_flux = -1
id_ustar_shelf = -1, id_shelf_mass = -1, id_mass_flux = -1, &
id_shelf_sfc_mass_flux = -1
!>@}

integer :: id_read_mass !< An integer handle used in time interpolation of
Expand Down Expand Up @@ -753,6 +754,8 @@ subroutine shelf_calc_flux(sfc_state_in, fluxes_in, Time, time_step, CS)
if (CS%id_shelf_mass > 0) call post_data(CS%id_shelf_mass, ISS%mass_shelf, CS%diag)
if (CS%id_area_shelf_h > 0) call post_data(CS%id_area_shelf_h, ISS%area_shelf_h, CS%diag)
if (CS%id_ustar_shelf > 0) call post_data(CS%id_ustar_shelf, fluxes%ustar_shelf, CS%diag)
if (CS%id_shelf_sfc_mass_flux > 0) call post_data(CS%id_shelf_sfc_mass_flux, fluxes%shelf_sfc_mass_flux, CS%diag)

if (CS%id_melt > 0) call post_data(CS%id_melt, fluxes%iceshelf_melt, CS%diag)
if (CS%id_thermal_driving > 0) call post_data(CS%id_thermal_driving, (sfc_state%sst-ISS%tfreeze), CS%diag)
if (CS%id_Sbdry > 0) call post_data(CS%id_Sbdry, Sbdry, CS%diag)
Expand Down Expand Up @@ -1815,6 +1818,8 @@ subroutine initialize_ice_shelf(param_file, ocn_grid, Time, CS, diag, forces_in,
if (CS%active_shelf_dynamics) then
CS%id_h_mask = register_diag_field('ice_shelf_model', 'h_mask', CS%diag%axesT1, CS%Time, &
'ice shelf thickness mask', 'none')
CS%id_shelf_sfc_mass_flux = register_diag_field('ice_shelf_model', 'sfc_mass_flux', CS%diag%axesT1, CS%Time, &
'ice shelf surface mass flux deposition from atmosphere', 'none')
endif
call MOM_IS_diag_mediator_close_registration(CS%diag)

Expand Down Expand Up @@ -1845,7 +1850,7 @@ subroutine initialize_ice_shelf_fluxes(CS, ocn_grid, US, fluxes_in)
! when SHELF_THERMO = True. These fluxes are necessary if one wants to
! use either ENERGETICS_SFC_PBL (ALE mode) or BULKMIXEDLAYER (layer mode).
call allocate_forcing_type(CS%Grid_in, fluxes_in, ustar=.true., shelf=.true., &
press=.true., water=CS%isthermo, heat=CS%isthermo)
press=.true., water=CS%isthermo, heat=CS%isthermo, shelf_sfc_accumulation = CS%active_shelf_dynamics)
else
call MOM_mesg("MOM_ice_shelf.F90, initialize_ice_shelf: allocating fluxes in solo mode.")
call allocate_forcing_type(CS%Grid_in, fluxes_in, ustar=.true., shelf=.true., press=.true.)
Expand Down

0 comments on commit c1a9131

Please sign in to comment.