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

Passing shelf_sfc_mass_flux to ice shelf (cleaned and rebased) #835

Closed
wants to merge 10 commits into from
9 changes: 9 additions & 0 deletions config_src/drivers/FMS_cap/MOM_surface_forcing_gfdl.F90
Original file line number Diff line number Diff line change
Expand Up @@ -201,6 +201,7 @@ module MOM_surface_forcing_gfdl
!! ice-shelves, expressed as a coefficient
!! for divergence damping, as determined
!! outside of the ocean model [m3 s-1]
real, pointer, dimension(:,:) :: shelf_sfc_mass_flux =>NULL() !< mass flux to surface of ice sheet [kg m-2 s-1]
integer :: xtype !< The type of the exchange - REGRID, REDIST or DIRECT
type(coupler_2d_bc_type) :: fluxes !< A structure that may contain an array of named fields
!! used for passive tracer fluxes.
Expand Down Expand Up @@ -464,6 +465,10 @@ subroutine convert_IOB_to_fluxes(IOB, fluxes, index_bounds, Time, valid_time, G,
call check_mask_val_consistency(IOB%calving(i-i0,j-j0), G%mask2dT(i,j), i, j, 'calving', G)
endif

if (associated(IOB%shelf_sfc_mass_flux)) then
fluxes%shelf_sfc_mass_flux(i,j) = kg_m2_s_conversion * IOB%shelf_sfc_mass_flux(i-i0,j-j0)
endif

if (associated(IOB%ustar_berg)) then
fluxes%ustar_berg(i,j) = US%m_to_Z*US%T_to_s * IOB%ustar_berg(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
Expand Down Expand Up @@ -1795,6 +1800,10 @@ subroutine ice_ocn_bnd_type_chksum(id, timestep, iobt)
chks = field_chksum( iobt%runoff ) ; if (root) write(outunit,100) 'iobt%runoff ', chks
chks = field_chksum( iobt%calving ) ; if (root) write(outunit,100) 'iobt%calving ', chks
chks = field_chksum( iobt%p ) ; if (root) write(outunit,100) 'iobt%p ', chks
if (associated(iobt%shelf_sfc_mass_flux)) then
chks = field_chksum( iobt%shelf_sfc_mass_flux ) ; if (root) write(outunit,100) 'iobt%shelf_sfc_mass_flux ',&
chks
endif
if (associated(iobt%ustar_berg)) then
chks = field_chksum( iobt%ustar_berg ) ; if (root) write(outunit,100) 'iobt%ustar_berg ', chks
endif
Expand Down
22 changes: 14 additions & 8 deletions src/ice_shelf/MOM_ice_shelf.F90
Original file line number Diff line number Diff line change
Expand Up @@ -856,7 +856,8 @@ subroutine shelf_calc_flux(sfc_state_in, fluxes_in, Time, time_step_in, CS)
if (CS%id_h_shelf > 0) call post_data(CS%id_h_shelf, ISS%h_shelf, CS%diag)
if (CS%id_dhdt_shelf > 0) call post_data(CS%id_dhdt_shelf, ISS%dhdt_shelf, CS%diag)
if (CS%id_h_mask > 0) call post_data(CS%id_h_mask,ISS%hmask,CS%diag)
call process_and_post_scalar_data(CS, vaf0, vaf0_A, vaf0_G, Itime_step, dh_adott, dh_bdott)
if (CS%active_shelf_dynamics) &
call process_and_post_scalar_data(CS, vaf0, vaf0_A, vaf0_G, Itime_step, dh_adott, dh_bdott)
call disable_averaging(CS%diag)

call cpu_clock_end(id_clock_shelf)
Expand Down Expand Up @@ -1875,7 +1876,12 @@ subroutine initialize_ice_shelf(param_file, ocn_grid, Time, CS, diag, Time_init,

CS%restart_output_dir = dirs%restart_output_dir


if (present(fluxes_in)) then
call initialize_ice_shelf_fluxes(CS, ocn_grid, US, fluxes_in)
call register_restart_field(fluxes_in%shelf_sfc_mass_flux, "sfc_mass_flux", .true., CS%restart_CSp, &
"ice shelf surface mass flux deposition from atmosphere", &
'kg m-2 s-1', conversion=US%RZ_T_to_kg_m2s)
endif

if (new_sim .and. (.not. (CS%override_shelf_movement .and. CS%mass_from_file))) then
! This model is initialized internally or from a file.
Expand Down Expand Up @@ -1998,11 +2004,12 @@ subroutine initialize_ice_shelf(param_file, ocn_grid, Time, CS, diag, Time_init,
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', conversion=1.0)
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', &
'kg m-2 s-1', conversion=US%RZ_T_to_kg_m2s)
endif

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', &
'kg m-2 s-1', conversion=US%RZ_T_to_kg_m2s)

! Scalars (area integrated over all ice sheets)
CS%id_vaf = register_scalar_field('ice_shelf_model', 'int_vaf', CS%diag%axesT1, CS%Time, &
'Area integrated ice sheet volume above floatation', 'm3', conversion=US%Z_to_m*US%L_to_m**2)
Expand Down Expand Up @@ -2178,7 +2185,6 @@ subroutine initialize_ice_shelf(param_file, ocn_grid, Time, CS, diag, Time_init,

call MOM_IS_diag_mediator_close_registration(CS%diag)

if (present(fluxes_in)) call initialize_ice_shelf_fluxes(CS, ocn_grid, US, fluxes_in)
if (present(forces_in)) call initialize_ice_shelf_forces(CS, ocn_grid, US, forces_in)

end subroutine initialize_ice_shelf
Expand Down Expand Up @@ -2351,8 +2357,8 @@ subroutine initialize_shelf_mass(G, param_file, CS, ISS, new_sim)

end subroutine initialize_shelf_mass
!> This subroutine applies net accumulation/ablation at the top surface to the dynamic ice shelf.
!>>acc_rate[m-s]=surf_mass_flux/density_ice is ablation/accumulation rate
!>>positive for accumulation negative for ablation
!! acc_rate[m-s]=surf_mass_flux/density_ice is ablation/accumulation rate
!! positive for accumulation negative for ablation
subroutine change_thickness_using_precip(CS, ISS, G, US, fluxes, time_step, Time)
type(ice_shelf_CS), intent(in) :: CS !< A pointer to the ice shelf control structure
type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure.
Expand Down