From 5d3d504ef3c9d4e6958aa6a8a32da66eea574f76 Mon Sep 17 00:00:00 2001 From: Robert Hallberg Date: Tue, 31 Dec 2024 19:24:37 -0500 Subject: [PATCH] (*)Conversion arguments for ice shelf scalar diags Revised volume_above_floatation and integrate_over_ice_sheet_area to return values in scaled units, and added conversion factors to the register_scalar_field calls in MOM_ice_shelf so that all unit unscaling of diagnostics occurs via the diag mediator and not in the code itself. After this commit, all of the dimensional scaling factors for diagnostics in the ice_shelf code occur via conversion factors that are immediately adjacent to the declaration of the units of those diagnostics, facilitating comparison for consistency. The declared units of one diagnostic, "taub_beta", were revised from "MPa s m-1" to "MPa yr m-1" to be consistent with its conversion factor, which was also corrected (essentially by a factor of [Z L-1 ~> 1]. Several missing unit conversion factors for rates of mass change were also added, and about 15 missing or incorrect units in lines documenting real variables were added or fixed. The input scale factor for the (perhaps unused) variable INPUT_VEL_ICE_SHELF was corrected from `US%m_s_to_L_T*US%m_to_Z` to just `US%m_s_to_L_T`. With these changes, all solutions are bitwise identical, apart from regional cases with a specified inflow when run with dimensional rescaling. The ice shelf diagnostics should now be invariant when dimensional rescaling is applied, and the units of the ice shelf diagnostics are now all consistent with the required rescaling factors. --- src/ice_shelf/MOM_ice_shelf.F90 | 187 ++++++++++++--------- src/ice_shelf/MOM_ice_shelf_dynamics.F90 | 27 +-- src/ice_shelf/MOM_ice_shelf_initialize.F90 | 4 +- 3 files changed, 124 insertions(+), 94 deletions(-) diff --git a/src/ice_shelf/MOM_ice_shelf.F90 b/src/ice_shelf/MOM_ice_shelf.F90 index 582518d4f1..ac12f18c9d 100644 --- a/src/ice_shelf/MOM_ice_shelf.F90 +++ b/src/ice_shelf/MOM_ice_shelf.F90 @@ -353,8 +353,8 @@ subroutine shelf_calc_flux(sfc_state_in, fluxes_in, Time, time_step_in, CS) character(len=160) :: mesg ! The text of an error message integer, dimension(2) :: EOSdom ! The i-computational domain for the equation of state integer :: i, j, is, ie, js, je, ied, jed, it1, it3 - real :: vaf0, vaf0_A, vaf0_G !The previous volumes above floatation [m3] - !for all ice sheets, Antarctica only, or Greenland only [m3] + real :: vaf0, vaf0_A, vaf0_G !The previous volumes above floatation [Z L2 ~> m3] + !for all ice sheets, Antarctica only, or Greenland only [Z L2 ~> m3] if (.not. associated(CS)) call MOM_error(FATAL, "shelf_calc_flux: "// & "initialize_ice_shelf must be called before shelf_calc_flux.") @@ -880,12 +880,12 @@ function integrate_over_ice_sheet_area(G, ISS, var, unscale, hemisphere) result( real, dimension(SZI_(G),SZJ_(G)), intent(in) :: var !< Ice variable to integrate in arbitrary units [A ~> a] real, intent(in) :: unscale !< Dimensional scaling for variable to integrate [a A-1 ~> 1] integer, optional, intent(in) :: hemisphere !< 0 for Antarctica only, 1 for Greenland only. Otherwise, all ice sheets - real :: var_out !< Variable integrated over the area of the ice sheet in arbitrary unscaled units [a m2] + real :: var_out !< Variable integrated over the area of the ice sheet in arbitrary scaled units [A L2 ~> a m2] ! Local variables integer :: IS_ID ! local copy of hemisphere real, dimension(SZI_(G),SZJ_(G)) :: var_cell !< Variable integrated over the ice-sheet area of each cell - !! in arbitrary units [a m2] + !! in arbitrary units [A L2 ~> a m2] integer, dimension(SZI_(G),SZJ_(G)) :: mask ! a mask for active cells depending on hemisphere indicated integer :: i,j @@ -913,7 +913,7 @@ function integrate_over_ice_sheet_area(G, ISS, var, unscale, hemisphere) result( if (mask(i,j)>0) var_cell(i,j) = var(i,j) * ISS%area_shelf_h(i,j) enddo; enddo - var_out = unscale*G%US%L_to_m**2 * reproducing_sum(var_cell, unscale=unscale*G%US%L_to_m**2) + var_out = reproducing_sum(var_cell, unscale=unscale*G%US%L_to_m**2) end function integrate_over_ice_sheet_area !> Converts the ice-shelf-to-ocean calving and calving_hflx variables from the ice-shelf state (ISS) type @@ -1853,7 +1853,7 @@ subroutine initialize_ice_shelf(param_file, ocn_grid, Time, CS, diag, Time_init, v_desc = var_desc("tauy_shelf", "Pa", "the meridional stress on the ocean under ice shelves", & hor_grid='Cv',z_grid='1') call register_restart_pair(sfc_state%taux_shelf, sfc_state%tauy_shelf, u_desc, v_desc, & - .false., CS%restart_CSp, conversion=US%RZ_T_to_kg_m2s*US%L_T_to_m_s) + .false., CS%restart_CSp, conversion=US%RLZ_T2_to_Pa) endif endif @@ -1997,134 +1997,161 @@ subroutine initialize_ice_shelf(param_file, ocn_grid, Time, CS, diag, Time_init, 'Fric vel under shelf', 'm/s', conversion=US%Z_to_m*US%s_to_T) 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') + '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 - !scalars (area integrated over all ice sheets) + ! 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') + 'Area integrated ice sheet volume above floatation', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_adott = register_scalar_field('ice_shelf_model', 'int_a', CS%diag%axesT1, CS%Time, & - 'Area integrated change in ice-sheet thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in ice-sheet thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_g_adott = register_scalar_field('ice_shelf_model', 'int_a_ground', CS%diag%axesT1, CS%Time, & - 'Area integrated change in grounded ice-sheet thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in grounded ice-sheet thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_f_adott = register_scalar_field('ice_shelf_model', 'int_a_float', CS%diag%axesT1, CS%Time, & - 'Area integrated change in floating ice-shelf thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in floating ice-shelf thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_bdott = register_scalar_field('ice_shelf_model', 'int_b', CS%diag%axesT1, CS%Time, & - 'Area integrated change in floating ice-shelf thickness '//& - 'due to basal accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in floating ice-shelf thickness '//& + 'due to basal accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_bdott_melt = register_scalar_field('ice_shelf_model', 'int_b_melt', CS%diag%axesT1, CS%Time, & - 'Area integrated basal melt over ice shelves during a DT_THERM time step', 'm3') + 'Area integrated basal melt over ice shelves during a DT_THERM time step', & + units='m3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_bdott_accum = register_scalar_field('ice_shelf_model', 'int_b_accum', CS%diag%axesT1, CS%Time, & - 'Area integrated basal accumulation over ice shelves during a DT_THERM a time step', 'm3') + 'Area integrated basal accumulation over ice shelves during a DT_THERM a time step', & + units='m3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_t_area = register_scalar_field('ice_shelf_model', 'tot_area', CS%diag%axesT1, CS%Time, & - 'Total ice-sheet area', 'm2') + 'Total ice-sheet area', 'm2', conversion=US%L_to_m**2) CS%id_f_area = register_scalar_field('ice_shelf_model', 'tot_area_float', CS%diag%axesT1, CS%Time, & - 'Total area of floating ice shelves', 'm2') + 'Total area of floating ice shelves', 'm2', conversion=US%L_to_m**2) CS%id_g_area = register_scalar_field('ice_shelf_model', 'tot_area_ground', CS%diag%axesT1, CS%Time, & - 'Total area of grounded ice sheets', 'm2') + 'Total area of grounded ice sheets', 'm2', conversion=US%L_to_m**2) !scalars (area integrated rates over all ice sheets) CS%id_dvafdt = register_scalar_field('ice_shelf_model', 'int_vafdot', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in ice-sheet volume above floatation', 'm3 s-1') - CS%id_adot = register_scalar_field('ice_shelf_model', 'int_adot', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in ice-sheet thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in ice-sheet volume above floatation', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) + CS%id_adot = register_scalar_field('ice_shelf_model', 'int_adot', CS%diag%axesT1, CS%Time, & + 'Area integrated rate of change in ice-sheet thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_g_adot = register_scalar_field('ice_shelf_model', 'int_adot_ground', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in grounded ice-sheet thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in grounded ice-sheet thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_f_adot = register_scalar_field('ice_shelf_model', 'int_adot_float', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in floating ice-shelf thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in floating ice-shelf thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_bdot = register_scalar_field('ice_shelf_model', 'int_bdot', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in ice-shelf thickness due to basal accum+melt', 'm3 s-1') + 'Area integrated rate of change in ice-shelf thickness due to basal accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_bdot_melt = register_scalar_field('ice_shelf_model', 'int_bdot_melt', CS%diag%axesT1, CS%Time, & - 'Area integrated basal melt rate over ice shelves', 'm3 s-1') + 'Area integrated basal melt rate over ice shelves', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_bdot_accum = register_scalar_field('ice_shelf_model', 'int_bdot_accum', CS%diag%axesT1, CS%Time, & - 'Area integrated basal accumulation rate over ice shelves', 'm3 s-1') + 'Area integrated basal accumulation rate over ice shelves', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) !scalars (area integrated over the Antarctic ice sheet) CS%id_Ant_vaf = register_scalar_field('ice_shelf_model', 'int_vaf_A', CS%diag%axesT1, CS%Time, & - 'Area integrated Antarctic ice sheet volume above floatation', 'm3') + 'Area integrated Antarctic ice sheet volume above floatation', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Ant_adott = register_scalar_field('ice_shelf_model', 'int_a_A', CS%diag%axesT1, CS%Time, & - 'Area integrated (Antarctic ice sheet) change in ice-sheet thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated (Antarctic ice sheet) change in ice-sheet thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Ant_g_adott = register_scalar_field('ice_shelf_model', 'int_a_ground_A', CS%diag%axesT1, CS%Time, & - 'Area integrated change in Antarctic grounded ice-sheet thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in Antarctic grounded ice-sheet thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Ant_f_adott = register_scalar_field('ice_shelf_model', 'int_a_float_A', CS%diag%axesT1, CS%Time, & - 'Area integrated change in Antarctic floating ice-shelf thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in Antarctic floating ice-shelf thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Ant_bdott = register_scalar_field('ice_shelf_model', 'int_b_A', CS%diag%axesT1, CS%Time, & - 'Area integrated change in Antarctic floating ice-shelf thickness '//& - 'due to basal accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in Antarctic floating ice-shelf thickness '//& + 'due to basal accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Ant_bdott_melt = register_scalar_field('ice_shelf_model', 'int_b_melt_A', CS%diag%axesT1, CS%Time, & - 'Area integrated basal melt over Antarctic ice shelves during a DT_THERM time step', 'm3') + 'Area integrated basal melt over Antarctic ice shelves during a DT_THERM time step', & + units='m3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Ant_bdott_accum = register_scalar_field('ice_shelf_model', 'int_b_accum_A', CS%diag%axesT1, CS%Time, & - 'Area integrated basal accumulation over Antarctic ice shelves during a DT_THERM a time step', 'm3') + 'Area integrated basal accumulation over Antarctic ice shelves during a DT_THERM a time step', & + units='m3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Ant_t_area = register_scalar_field('ice_shelf_model', 'tot_area_A', CS%diag%axesT1, CS%Time, & - 'Total area of Antarctic ice sheet', 'm2') + 'Total area of Antarctic ice sheet', 'm2', conversion=US%L_to_m**2) CS%id_Ant_f_area = register_scalar_field('ice_shelf_model', 'tot_area_float_A', CS%diag%axesT1, CS%Time, & - 'Total area of Antarctic floating ice shelves', 'm2') + 'Total area of Antarctic floating ice shelves', 'm2', conversion=US%L_to_m**2) CS%id_Ant_g_area = register_scalar_field('ice_shelf_model', 'tot_area_ground_A', CS%diag%axesT1, CS%Time, & - 'Total area of Antarctic grounded ice sheet', 'm2') + 'Total area of Antarctic grounded ice sheet', 'm2', conversion=US%L_to_m**2) !scalars (area integrated rates over the Antarctic ice sheet) CS%id_Ant_dvafdt = register_scalar_field('ice_shelf_model', 'int_vafdot_A', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Antarctic ice-sheet volume above floatation', 'm3 s-1') - CS%id_Ant_adot = register_scalar_field('ice_shelf_model', 'int_adot_A', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Antarctic ice-sheet thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in Antarctic ice-sheet volume above floatation', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) + CS%id_Ant_adot = register_scalar_field('ice_shelf_model', 'int_adot_A', CS%diag%axesT1, CS%Time, & + 'Area integrated rate of change in Antarctic ice-sheet thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Ant_g_adot = register_scalar_field('ice_shelf_model', 'int_adot_ground_A', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Antarctic grounded ice-sheet thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in Antarctic grounded ice-sheet thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Ant_f_adot = register_scalar_field('ice_shelf_model', 'int_adot_float_A', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Antarctic floating ice-shelf thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in Antarctic floating ice-shelf thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Ant_bdot = register_scalar_field('ice_shelf_model', 'int_bdot_A', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Antarctic ice-shelf thickness due to basal accum+melt', 'm3 s-1') + 'Area integrated rate of change in Antarctic ice-shelf thickness due to basal accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Ant_bdot_melt = register_scalar_field('ice_shelf_model', 'int_bdot_melt_A', CS%diag%axesT1, CS%Time, & - 'Area integrated basal melt rate over Antarctic ice shelves', 'm3 s-1') + 'Area integrated basal melt rate over Antarctic ice shelves', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Ant_bdot_accum = register_scalar_field('ice_shelf_model', 'int_bdot_accum_A', CS%diag%axesT1, CS%Time, & - 'Area integrated basal accumulation rate over Antarctic ice shelves', 'm3 s-1') + 'Area integrated basal accumulation rate over Antarctic ice shelves', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) !scalars (area integrated over the Greenland ice sheet) CS%id_Gr_vaf = register_scalar_field('ice_shelf_model', 'int_vaf_G', CS%diag%axesT1, CS%Time, & - 'Area integrated Greenland ice sheet volume above floatation', 'm3') + 'Area integrated Greenland ice sheet volume above floatation', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Gr_adott = register_scalar_field('ice_shelf_model', 'int_a_G', CS%diag%axesT1, CS%Time, & - 'Area integrated (Greenland ice sheet) change in ice-sheet thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated (Greenland ice sheet) change in ice-sheet thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Gr_g_adott = register_scalar_field('ice_shelf_model', 'int_a_ground_G', CS%diag%axesT1, CS%Time, & - 'Area integrated change in Greenland grounded ice-sheet thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in Greenland grounded ice-sheet thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Gr_f_adott = register_scalar_field('ice_shelf_model', 'int_a_float_G', CS%diag%axesT1, CS%Time, & - 'Area integrated change in Greenland floating ice-shelf thickness ' //& - 'due to surface accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in Greenland floating ice-shelf thickness ' //& + 'due to surface accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Gr_bdott = register_scalar_field('ice_shelf_model', 'int_b_G', CS%diag%axesT1, CS%Time, & - 'Area integrated change in Greenland floating ice-shelf thickness '//& - 'due to basal accum+melt during a DT_THERM time step', 'm3') + 'Area integrated change in Greenland floating ice-shelf thickness '//& + 'due to basal accum+melt during a DT_THERM time step', 'm3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Gr_bdott_melt = register_scalar_field('ice_shelf_model', 'int_b_melt_G', CS%diag%axesT1, CS%Time, & - 'Area integrated basal melt over Greenland ice shelves during a DT_THERM time step', 'm3') + 'Area integrated basal melt over Greenland ice shelves during a DT_THERM time step', & + units='m3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Gr_bdott_accum = register_scalar_field('ice_shelf_model', 'int_b_accum_G', CS%diag%axesT1, CS%Time, & - 'Area integrated basal accumulation over Greenland ice shelves during a DT_THERM a time step', 'm3') + 'Area integrated basal accumulation over Greenland ice shelves during a DT_THERM a time step', & + units='m3', conversion=US%Z_to_m*US%L_to_m**2) CS%id_Gr_t_area = register_scalar_field('ice_shelf_model', 'tot_area_G', CS%diag%axesT1, CS%Time, & - 'Total area of Greenland ice sheet', 'm2') + 'Total area of Greenland ice sheet', 'm2', conversion=US%L_to_m**2) CS%id_Gr_f_area = register_scalar_field('ice_shelf_model', 'tot_area_float_G', CS%diag%axesT1, CS%Time, & - 'Total area of Greenland floating ice shelves', 'm2') + 'Total area of Greenland floating ice shelves', 'm2', conversion=US%L_to_m**2) CS%id_Gr_g_area = register_scalar_field('ice_shelf_model', 'tot_area_ground_G', CS%diag%axesT1, CS%Time, & - 'Total area of Greenland grounded ice sheet', 'm2') + 'Total area of Greenland grounded ice sheet', 'm2', conversion=US%L_to_m**2) !scalars (area integrated rates over the Greenland ice sheet) CS%id_Gr_dvafdt = register_scalar_field('ice_shelf_model', 'int_vafdot_G', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Greenland ice-sheet volume above floatation', 'm3 s-1') - CS%id_Gr_adot = register_scalar_field('ice_shelf_model', 'int_adot_G', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Greenland ice-sheet thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in Greenland ice-sheet volume above floatation', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) + CS%id_Gr_adot = register_scalar_field('ice_shelf_model', 'int_adot_G', CS%diag%axesT1, CS%Time, & + 'Area integrated rate of change in Greenland ice-sheet thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Gr_g_adot = register_scalar_field('ice_shelf_model', 'int_adot_ground_G', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Greenland grounded ice-sheet thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in Greenland grounded ice-sheet thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Gr_f_adot = register_scalar_field('ice_shelf_model', 'int_adot_float_G', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Greenland floating ice-shelf thickness due to surface accum+melt', 'm3 s-1') + 'Area integrated rate of change in Greenland floating ice-shelf thickness due to surface accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Gr_bdot = register_scalar_field('ice_shelf_model', 'int_bdot_G', CS%diag%axesT1, CS%Time, & - 'Area integrated rate of change in Greenland ice-shelf thickness due to basal accum+melt', 'm3 s-1') + 'Area integrated rate of change in Greenland ice-shelf thickness due to basal accum+melt', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Gr_bdot_melt = register_scalar_field('ice_shelf_model', 'int_bdot_melt_G', CS%diag%axesT1, CS%Time, & - 'Area integrated basal melt rate over Greenland ice shelves', 'm3 s-1') + 'Area integrated basal melt rate over Greenland ice shelves', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) CS%id_Gr_bdot_accum = register_scalar_field('ice_shelf_model', 'int_bdot_accum_G', CS%diag%axesT1, CS%Time, & - 'Area integrated basal accumulation rate over Greenland ice shelves', 'm3 s-1') + 'Area integrated basal accumulation rate over Greenland ice shelves', & + units='m3 s-1', conversion=US%Z_to_m*US%L_to_m**2*US%s_to_T) !Flags to calculate diagnostics related to surface/basal mass balance if (CS%id_adott>0 .or. CS%id_g_adott>0 .or. CS%id_f_adott>0 .or. & @@ -2520,7 +2547,7 @@ subroutine solo_step_ice_shelf(CS, time_interval, nsteps, Time, min_time_step_in ! coupled ice-ocean dynamics. integer :: is, ie, js, je, i, j real :: vaf0, vaf0_A, vaf0_G !The previous volumes above floatation - !for all ice sheets, Antarctica only, or Greenland only [m3] + !for all ice sheets, Antarctica only, or Greenland only [Z L2 ~> m3] real, dimension(SZI_(CS%grid),SZJ_(CS%grid)) :: & dh_adott_sum, & ! Surface melt/accumulation over a full time step, used for diagnostics [Z ~> m] dh_adott ! Surface melt/accumulation over a partial time step, used for diagnostics [Z ~> m] @@ -2604,17 +2631,19 @@ end subroutine solo_step_ice_shelf !> Post_data calls for ice-sheet scalars subroutine process_and_post_scalar_data(CS, vaf0, vaf0_A, vaf0_G, Itime_step, dh_adott, dh_bdott) type(ice_shelf_CS), pointer :: CS !< A pointer to the ice shelf control structure - real :: vaf0 !< The previous volumes above floatation for all ice sheets [m3] - real :: vaf0_A !< The previous volumes above floatation for the Antarctic ice sheet [m3] - real :: vaf0_G !< The previous volumes above floatation for the Greenland ice sheet [m3] + real :: vaf0 !< The previous volumes above floatation for all ice sheets [Z L2 ~> m3] + real :: vaf0_A !< The previous volumes above floatation for the Antarctic ice sheet [Z L2 ~> m3] + real :: vaf0_G !< The previous volumes above floatation for the Greenland ice sheet [Z L2 ~> m3] real :: Itime_step !< Inverse of the time step [T-1 ~> s-1] real, dimension(SZI_(CS%grid),SZJ_(CS%grid)) :: dh_adott !< Surface (plus basal if solo shelf mode) !! melt/accumulation over a time step [Z ~> m] real, dimension(SZI_(CS%grid),SZJ_(CS%grid)) :: dh_bdott !< Surface (plus basal if solo shelf mode) !! melt/accumulation over a time step [Z ~> m] + + ! Local variables real, dimension(SZI_(CS%grid),SZJ_(CS%grid)) :: tmp ! Temporary field used when calculating diagnostics [various] real, dimension(SZI_(CS%grid),SZJ_(CS%grid)) :: ones ! Temporary field used when calculating diagnostics [various] - real :: vaf ! The current ice-sheet volume above floatation [m3] + real :: vaf ! The current ice-sheet volume above floatation [Z L2 ~> m3] real :: val ! Temporary value when calculating scalar diagnostics [various] type(ocean_grid_type), pointer :: G => NULL() ! A pointer to the ocean's grid structure type(unit_scale_type), pointer :: US => NULL() ! Pointer to a structure containing various unit conversion factors diff --git a/src/ice_shelf/MOM_ice_shelf_dynamics.F90 b/src/ice_shelf/MOM_ice_shelf_dynamics.F90 index 140b9746d8..09eea05127 100644 --- a/src/ice_shelf/MOM_ice_shelf_dynamics.F90 +++ b/src/ice_shelf/MOM_ice_shelf_dynamics.F90 @@ -108,10 +108,10 @@ module MOM_ice_shelf_dynamics real, pointer, dimension(:,:) :: basal_traction => NULL() !< The area-integrated taub_beta field !! (m2 Pa s m-1, or kg s-1) related to the nonlinear part - !! of "linearized" basal stress (Pa) [R L3 T-1 ~> kg s-1] + !! of "linearized" basal stress (Pa) [R Z L2 T-1 ~> kg s-1] !! The exact form depends on basal law exponent and/or whether flow is "hybridized" a la Goldberg 2011 real, pointer, dimension(:,:) :: C_basal_friction => NULL()!< Coefficient in sliding law tau_b = C u^(n_basal_fric), - !! units= Pa (s m-1)^(n_basal_fric) + !! units= [Pa (s m-1)^(n_basal_fric)] real, pointer, dimension(:,:) :: OD_rt => NULL() !< A running total for calculating OD_av [Z ~> m]. real, pointer, dimension(:,:) :: ground_frac_rt => NULL() !< A running total for calculating ground_frac. real, pointer, dimension(:,:) :: OD_av => NULL() !< The time average open ocean depth [Z ~> m]. @@ -169,7 +169,7 @@ module MOM_ice_shelf_dynamics !! i.e. dt <= CFL_factor * min(dx / u) [nondim] real :: min_h_shelf !< The minimum ice thickness used during ice dynamics [L ~> m]. - real :: min_basal_traction !< The minimum basal traction for grounded ice (Pa m-1 s) [R L T-1 ~> kg m-2 s-1] + real :: min_basal_traction !< The minimum basal traction for grounded ice (Pa m-1 s) [R Z T-1 ~> kg m-2 s-1] real :: max_surface_slope !< The maximum allowed ice-sheet surface slope (to ignore, set to zero) [nondim] real :: min_ice_visc !< The minimum allowed Glen's law ice viscosity (Pa s), in [R L2 T-1 ~> kg m-1 s-1]. @@ -358,7 +358,7 @@ subroutine register_ice_shelf_dyn_restarts(G, US, param_file, CS, restart_CS) allocate(CS%t_shelf(isd:ied,jsd:jed), source=T_shelf_missing) ! [C ~> degC] allocate(CS%ice_visc(isd:ied,jsd:jed,CS%visc_qps), source=0.0) allocate(CS%AGlen_visc(isd:ied,jsd:jed), source=2.261e-25) ! [Pa-3 s-1] - allocate(CS%basal_traction(isd:ied,jsd:jed), source=0.0) ! [R L3 T-1 ~> kg s-1] + allocate(CS%basal_traction(isd:ied,jsd:jed), source=0.0) ! [R Z L2 T-1 ~> kg s-1] allocate(CS%C_basal_friction(isd:ied,jsd:jed), source=5.0e10) ! [Pa (s m-1)^n_sliding] allocate(CS%OD_av(isd:ied,jsd:jed), source=0.0) allocate(CS%ground_frac(isd:ied,jsd:jed), source=0.0) @@ -838,7 +838,7 @@ subroutine initialize_ice_shelf_dyn(param_file, Time, ISS, CS, G, US, diag, new_ CS%id_visc_shelf = register_diag_field('ice_shelf_model','ice_visc',CS%diag%axesT1, Time, & 'vi-viscosity', 'Pa m s', conversion=US%RL2_T2_to_Pa*US%Z_to_m*US%T_to_s) !vertically integrated viscosity CS%id_taub = register_diag_field('ice_shelf_model','taub_beta',CS%diag%axesT1, Time, & - 'taub', 'MPa s m-1', conversion=1e-6*US%RL2_T2_to_Pa/(365.0*86400.0*US%L_T_to_m_s)) + 'taub', units='MPa yr m-1', conversion=1e-6*US%RLZ_T2_to_Pa/(365.0*86400.0*US%L_T_to_m_s)) CS%id_OD_av = register_diag_field('ice_shelf_model','OD_av',CS%diag%axesT1, Time, & 'intermediate ocean column thickness passed to ice model', 'm', conversion=US%Z_to_m) @@ -1010,10 +1010,10 @@ subroutine volume_above_floatation(CS, G, ISS, vaf, hemisphere) type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf. type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe !! the ice-shelf state - real, intent(out) :: vaf !< area integrated volume above floatation [m3] + real, intent(out) :: vaf !< area integrated volume above floatation [Z L2 ~> m3] integer, optional, intent(in) :: hemisphere !< 0 for Antarctica only, 1 for Greenland only. Otherwise, all ice sheets integer :: IS_ID ! local copy of hemisphere - real, dimension(SZI_(G),SZJ_(G)) :: vaf_cell !< cell-wise volume above floatation [m3] + real, dimension(SZI_(G),SZJ_(G)) :: vaf_cell !< cell-wise volume above floatation [Z L2 ~> m3] integer, dimension(SZI_(G),SZJ_(G)) :: mask ! a mask for active cells depending on hemisphere indicated integer :: is,ie,js,je,i,j real :: rhoi_rhow, rhow_rhoi @@ -1057,7 +1057,7 @@ subroutine volume_above_floatation(CS, G, ISS, vaf, hemisphere) endif enddo; enddo - vaf = G%US%Z_to_m*G%US%L_to_m**2 * reproducing_sum(vaf_cell, unscale=G%US%Z_to_m*G%US%L_to_m**2) + vaf = reproducing_sum(vaf_cell, unscale=G%US%Z_to_m*G%US%L_to_m**2) end subroutine volume_above_floatation !> multiplies a variable with the ice sheet grounding fraction @@ -2642,7 +2642,7 @@ subroutine CG_action(CS, uret, vret, u_shlf, v_shlf, Phi, Phisub, umask, vmask, !! relative to sea-level [Z ~> m]. real, dimension(SZDI_(G),SZDJ_(G)), & intent(in) :: basal_trac !< Area-integrated taub_beta field related to the nonlinear - !! part of the "linearized" basal stress [R L3 T-1 ~> kg s-1]. + !! part of the "linearized" basal stress [R Z L2 T-1 ~> kg s-1]. real, intent(in) :: dens_ratio !< The density of ice divided by the density !! of seawater, nondimensional @@ -2675,10 +2675,11 @@ subroutine CG_action(CS, uret, vret, u_shlf, v_shlf, Phi, Phisub, umask, vmask, real :: uq, vq ! Interpolated velocities [L T-1 ~> m s-1] integer :: iq, jq, iphi, jphi, i, j, ilq, jlq, Itgt, Jtgt, qp, qpv logical :: visc_qp4 - real, dimension(2) :: xquad - real, dimension(2,2) :: Ucell, Vcell, Hcell, Usub, Vsub - real, dimension(2,2,4) :: uret_qp, vret_qp - real, dimension(SZDIB_(G),SZDJB_(G),4) :: uret_b, vret_b + real, dimension(2) :: xquad ! Nondimensional quadrature ratios [nondim] + real, dimension(2,2) :: Ucell, Vcell, Usub, Vsub ! Velocities at the nodal points around the cell [L T-1 ~> m s-1] + real, dimension(2,2) :: Hcell ! Ice shelf thickness at notal (corner) points [Z ~> m] + real, dimension(2,2,4) :: uret_qp, vret_qp ! Temporary arrays in [R Z L3 T-2 ~> kg m s-2] + real, dimension(SZDIB_(G),SZDJB_(G),4) :: uret_b, vret_b ! Temporary arrays in [R Z L3 T-2 ~> kg m s-2] xquad(1) = .5 * (1-sqrt(1./3)) ; xquad(2) = .5 * (1+sqrt(1./3)) diff --git a/src/ice_shelf/MOM_ice_shelf_initialize.F90 b/src/ice_shelf/MOM_ice_shelf_initialize.F90 index f976187c2b..3dfe1045cf 100644 --- a/src/ice_shelf/MOM_ice_shelf_initialize.F90 +++ b/src/ice_shelf/MOM_ice_shelf_initialize.F90 @@ -321,7 +321,7 @@ subroutine initialize_ice_shelf_boundary_channel(u_face_mask_bdry, v_face_mask_b call get_param(PF, mdl, "INPUT_VEL_ICE_SHELF", input_vel, & "inflow ice velocity at upstream boundary", & - units="m s-1", default=0., scale=US%m_s_to_L_T*US%m_to_Z) !### This conversion factor is wrong? + units="m s-1", default=0., scale=US%m_s_to_L_T) call get_param(PF, mdl, "INPUT_THICK_ICE_SHELF", input_thick, & "flux thickness at upstream boundary", & units="m", default=1000., scale=US%m_to_Z) @@ -556,7 +556,7 @@ end subroutine initialize_ice_shelf_boundary_from_file subroutine initialize_ice_C_basal_friction(C_basal_friction, G, US, PF) type(ocean_grid_type), intent(in) :: G !< The ocean's grid structure real, dimension(SZDI_(G),SZDJ_(G)), & - intent(inout) :: C_basal_friction !< Ice-stream basal friction + intent(inout) :: C_basal_friction !< Ice-stream basal friction [Pa (s m-1)^(n_basal_fric)] type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors type(param_file_type), intent(in) :: PF !< A structure to parse for run-time parameters