mom_ice_shelf_dynamics module reference¶
Implements a crude placeholder for a later implementation of full ice shelf dynamics.
Data Types¶
The control structure for the ice shelf dynamics. |
|
A container for loop bounds. |
Functions/Subroutines¶
used for flux limiting in advective subroutines Van Leer limiter (source: Wikipedia) The return value is between 0 and 2 [nondim]. |
|
Calculate area of quadrilateral. |
|
This subroutine is used to register any fields related to the ice shelf dynamics that should be written to or read from the restart file. |
|
Initializes shelf model data, parameters and diagnostics. |
|
This function returns the global maximum advective timestep that can be taken based on the current ice velocities. |
|
This subroutine updates the ice shelf velocities, mass, stresses and properties due to the ice shelf dynamics. |
|
multiplies a variable with the ice sheet grounding fraction |
|
Ice shelf dynamics post_data calls. |
|
Calculate cell-centered, area-averaged, vertically integrated ice viscosity for diagnostics. |
|
Writes the total ice shelf kinetic energy and mass to an ascii file. |
|
This subroutine takes the velocity (on the Bgrid) and timesteps h_t = - div (uh) once. |
|
This subroutine computes u- and v-velocities of the ice shelf iterating on non-linear ice viscosity. |
|
Apply a very simple calving law using a minimum thickness rule. |
|
Calculate driving stress using cell-centered bed elevation and ice thickness. |
|
returns the diagonal entries of the matrix for a Jacobi preconditioning |
|
Post_data calls related to ice-sheet flux divergence, strain-rate, and deviatoric stress. |
|
Update depth integrated viscosity, based on horizontal strain rates. |
|
Update basal shear. |
|
This subroutine calculates the gradients of bilinear basis elements that that are centered at the vertices of the cell. |
|
This subroutine calculates the gradients of bilinear basis elements that are centered at the vertices of the cell using a locally orthogoal MOM6 grid. |
|
This subroutine calculates the gradients of bilinear basis elements that are centered at the vertices of the cell using a locally orthogoal MOM6 grid. |
|
Interpolate the ice shelf thickness from tracer point to nodal points, subject to a mask. |
|
Deallocates all memory associated with the ice shelf dynamics module. |
|
This subroutine updates the vertically averaged ice shelf temperature. |
|
Detailed Description¶
Implements a crude placeholder for a later implementation of full ice shelf dynamics.
Type Documentation¶
-
type
mom_ice_shelf_dynamics/ice_shelf_dyn_cs¶ The control structure for the ice shelf dynamics.
- Type fields:
%id_u_shelf[integer] :: Diagnostic handles.%id_v_shelf[integer] :: Diagnostic handles.%id_shelf_speed[integer] :: Diagnostic handles.%id_t_shelf[integer] :: Diagnostic handles.%id_taudx_shelf[integer] :: Diagnostic handles.%id_taudy_shelf[integer] :: Diagnostic handles.%id_taud_shelf[integer] :: Diagnostic handles.%id_bed_elev[integer] :: Diagnostic handles.%id_ground_frac[integer] :: Diagnostic handles.%id_col_thick[integer] :: Diagnostic handles.%id_od_av[integer] :: Diagnostic handles.%id_float_cond[integer] :: Diagnostic handles.%id_u_mask[integer] :: Diagnostic handles.%id_v_mask[integer] :: Diagnostic handles.%id_ufb_mask[integer] :: Diagnostic handles.%id_vfb_mask[integer] :: Diagnostic handles.%id_t_mask[integer] :: Diagnostic handles.%id_sx_shelf[integer] :: Diagnostic handles.%id_sy_shelf[integer] :: Diagnostic handles.%id_surf_slope_mag_shelf[integer] :: Diagnostic handles.%id_duhdx[integer] :: Diagnostic handles.%id_dvhdy[integer] :: Diagnostic handles.%id_fluxdiv[integer] :: Diagnostic handles.%id_strainrate_xx[integer] :: Diagnostic handles.%id_strainrate_yy[integer] :: Diagnostic handles.%id_strainrate_xy[integer] :: Diagnostic handles.%id_pstrainrate_1[integer] :: Diagnostic handles.%id_pstrainrate_2[integer] :: Diagnostic handles.%id_devstress_xx[integer] :: Diagnostic handles.%id_devstress_yy[integer] :: Diagnostic handles.%id_devstress_xy[integer] :: Diagnostic handles.%id_pdevstress_1[integer] :: Diagnostic handles.%id_pdevstress_2[integer] :: Diagnostic handles.%id_h_after_uflux[integer] :: Diagnostic handles for debugging.%id_h_after_vflux[integer] :: Diagnostic handles for debugging.%id_h_after_adv[integer] :: Diagnostic handles for debugging.%id_visc_shelf[integer] :: Diagnostic handles for debugging.%id_taub[integer] :: Diagnostic handles for debugging.%u_shelf[real(:,:),pointer] :: the zonal velocity of the ice shelf/sheet on q-points (B grid) [L T-1 ~> m s-1]%v_shelf[real(:,:),pointer] :: the meridional velocity of the ice shelf/sheet on q-points (B grid) [L T-1 ~> m s-1]%taudx_shelf[real(:,:),pointer] :: the zonal driving stress of the ice shelf/sheet on q-points (C grid) [R L2 T-2 ~> Pa]%taudy_shelf[real(:,:),pointer] :: the meridional driving stress of the ice shelf/sheet on q-points (C grid) [R L2 T-2 ~> Pa]%sx_shelf[real(:,:),pointer] :: the zonal surface slope of the ice shelf/sheet on q-points (B grid) [nondim]%sy_shelf[real(:,:),pointer] :: the meridional surface slope of the ice shelf/sheet on q-points (B grid) [nondim]%u_face_mask[real(:,:),pointer] :: mask for velocity boundary conditions on the C-grid u-face - this is because the FEM cares about FACES THAT GET INTEGRATED OVER, not vertices. Will represent boundary conditions on computational boundary (or permanent boundary between fast-moving and near-stagnant ice FOR NOW: 1=interior bdry, 0=no-flow boundary, 2=stress bdry condition, 3=inhomogeneous Dirichlet boundary for u and v, 4=flux boundary: at these faces a flux will be specified which will override velocities; a homogeneous velocity condition will be specified (this seems to give the solver less difficulty) 5=inhomogenous Dirichlet boundary for u only. 6=inhomogenous Dirichlet boundary for v only%v_face_mask[real(:,:),pointer] :: A mask for velocity boundary conditions on the C-grid v-face, with valued defined similarly to u_face_mask, but 5 is Dirichlet for v and 6 is Dirichlet for u.%u_face_mask_bdry[real(:,:),pointer] :: A duplicate copy of u_face_mask?%v_face_mask_bdry[real(:,:),pointer] :: A duplicate copy of v_face_mask?%u_flux_bdry_val[real(:,:),pointer] :: The ice volume flux per unit face length into the cell through open boundary u-faces (where u_face_mask=4) [Z L T-1 ~> m2 s-1].%v_flux_bdry_val[real(:,:),pointer] :: The ice volume flux per unit face length into the cell through open boundary v-faces (where v_face_mask=4) [Z L T-1 ~> m2 s-1]??%umask[real(:,:),pointer] :: u-mask on the actual degrees of freedom (B grid) 1=normal node, 3=inhomogeneous boundary node, 0 - no flow node (will also get ice-free nodes)%vmask[real(:,:),pointer] :: v-mask on the actual degrees of freedom (B grid) 1=normal node, 3=inhomogeneous boundary node, 0 - no flow node (will also get ice-free nodes)%calve_mask[real(:,:),pointer] :: a mask to prevent the ice shelf front from advancing past its initial position (but it may retreat)%t_shelf[real(:,:),pointer] :: Vertically integrated temperature in the ice shelf/stream, on corner-points (B grid) [C ~> degC].%tmask[real(:,:),pointer] :: A mask on tracer points that is 1 where there is ice.%ice_visc[real(:,:,:),pointer] :: Area and depth-integrated Glen’s law ice viscosity (Pa m3 s) in [R L4 Z T-1 ~> kg m2 s-1]. at either 1 (cell-centered) or 4 quadrature points per cell.%aglen_visc[real(:,:),pointer] :: Ice-stiffness parameter in Glen’s law ice viscosity, often in [Pa-3 s-1] if n_Glen is 3.%u_bdry_val[real(:,:),pointer] :: The zonal ice velocity at inflowing boundaries [L yr-1 ~> m yr-1].%v_bdry_val[real(:,:),pointer] :: The meridional ice velocity at inflowing boundaries [L yr-1 ~> m yr-1].%h_bdry_val[real(:,:),pointer] :: The ice thickness at inflowing boundaries [Z ~> m].%t_bdry_val[real(:,:),pointer] :: The ice temperature at inflowing boundaries [C ~> degC].%bed_elev[real(:,:),pointer] :: The bed elevation used for ice dynamics [Z ~> m], relative to mean sea-level. This is the same as GbathyT+Z_ref, when below sea-level. Sign convention: positive below sea-level, negative above.%basal_traction[real(:,:),pointer] :: 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 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.%c_basal_friction[real(:,:),pointer] :: Coefficient in sliding law tau_b = C u^(n_basal_fric), units of [R L Z T-2 (s m-1)^(n_basal_fric) ~> Pa (s m-1)^(n_basal_fric)].%od_rt[real(:,:),pointer] :: A running total for calculating OD_av [Z ~> m].%ground_frac_rt[real(:,:),pointer] :: A running total for calculating ground_frac.%od_av[real(:,:),pointer] :: The time average open ocean depth [Z ~> m].%ground_frac[real(:,:),pointer] :: Fraction of the time a cell is “exposed”, i.e. the column thickness is below a threshold and interacting with the rock [nondim]. When this is 1, the ice-shelf is grounded.%float_cond[real(:,:),pointer] :: If GL_regularize=true, indicates cells containing the grounding line (float_cond=1) or not (float_cond=0)%phi[real(:,:,:,:),pointer] :: The gradients of bilinear basis elements at Gaussian 4 quadrature points surrounding the cell vertices [L-1 ~> m-1].%phic[real(:,:,:),pointer] :: The gradients of bilinear basis elements at 1 cell-centered quadrature point per cell [L-1 ~> m-1].%phisub[real(:,:,:,:,:,:),pointer] :: Quadrature structure weights at subgridscale locations for finite element calculations [nondim].%od_rt_counter[integer] :: A counter of the number of contributions to OD_rt.%velocity_update_time_step[real] :: The time interval over which to update the ice shelf velocity using the nonlinear elliptic equation, or 0 to update every timestep [T ~> s].%elapsed_velocity_time[real] :: The elapsed time since the ice velocities were last updated [T ~> s].%g_earth[real] :: The gravitational acceleration [L2 Z-1 T-2 ~> m s-2].%density_ice[real] :: A typical density of ice [R ~> kg m-3].%cp_ice[real] :: The heat capacity of fresh ice [Q C-1 ~> J kg-1 degC-1].%advect_shelf[logical] :: If true (default), advect ice shelf and evolve thickness.%reentrant_x[logical] :: If true, the domain is zonally reentrant.%reentrant_y[logical] :: If true, the domain is meridionally reentrant.%alternate_first_direction_is[logical] :: If true, alternate whether the x- or y-direction updates occur first in directionally split parts of the calculation.%first_direction_is[integer] :: An integer that indicates which direction is to be updated first in directionally split parts of the ice sheet calculation (e.g. advection).%first_dir_restart_is[real] :: A real copy of CSfirst_direction_IS for use in restart files.%visc_qps[integer] :: The number of quadrature points per cell (1 or 4) on which to calculate ice viscosity.%ice_viscosity_compute[character (len=40)] :: Specifies whether the ice viscosity is computed internally according to Glen’s flow law; is constant (for debugging purposes) or using observed strain rates and read from a file.%gl_regularize[logical] :: Specifies whether to regularize the floatation condition at the grounding line as in Goldberg Holland Schoof 2009.%n_sub_regularize[integer] :: partition of cell over which to integrate for interpolated grounding line the (rectangular) is divided into nxn equally-sized rectangles, over which basal contribution is integrated (iterative quadrature)%gl_couple[logical] :: whether to let the floatation condition be determined by ocean column thickness means update_OD_ffrac will be called (note: GL_regularize and GL_couple should be exclusive)%cfl_factor[real] :: A factor used to limit subcycled advective timestep in uncoupled runs i.e. dt <= CFL_factor * min(dx / u) [nondim].%min_h_shelf[real] :: The minimum ice thickness used during ice dynamics [Z ~> m].%min_basal_traction[real] :: The minimum basal traction for grounded ice (Pa m-1 s) [R Z T-1 ~> kg m-2 s-1].%max_surface_slope[real] :: The maximum allowed ice-sheet surface slope (to ignore, set to zero) [nondim].%min_ice_visc[real] :: The minimum allowed Glen’s law ice viscosity (Pa s), in [R L2 T-1 ~> kg m-1 s-1].%n_glen[real] :: Nonlinearity exponent in Glen’s Law [nondim].%eps_glen_min[real] :: Min. strain rate to avoid infinite Glen’s law viscosity, [T-1 ~> s-1].%n_basal_fric[real] :: Exponent in sliding law tau_b = C u^(m_slide) [nondim].%coulombfriction[logical] :: Use Coulomb friction law (Schoof 2005, Gagliardini et al 2007)%cf_minn[real] :: Minimum Coulomb friction effective pressure [R Z L T-2 ~> Pa].%cf_postpeak[real] :: Coulomb friction post peak exponent [nondim].%cf_max[real] :: Coulomb friction maximum coefficient [nondim].%density_ocean_avg[real] :: A typical ocean density [R ~> kg m-3]. This does not affect ocean circulation or thermodynamics. It is used to estimate the gravitational driving force at the shelf front (until we think of a better way to do it, but any difference will be negligible).%thresh_float_col_depth[real] :: The water column depth over which the shelf if considered to be floating.%moving_shelf_front[logical] :: Specify whether to advance shelf front (and calve).%calve_to_mask[logical] :: If true, calve off the ice shelf when it passes the edge of a mask.%min_thickness_simple_calve[real] :: min. ice shelf thickness criteria for calving [Z ~> m].%t_shelf_missing[real] :: An ice shelf temperature to use where there is no ice shelf [C ~> degC].%cg_tolerance[real] :: The tolerance in the CG solver, relative to initial residual, that determines when to stop the conjugate gradient iterations [nondim].%nonlinear_tolerance[real] :: The fractional nonlinear tolerance, relative to the initial error, that sets when to stop the iterative velocity solver [nondim].%cg_max_iterations[integer] :: The maximum number of iterations that can be used in the CG solver.%nonlin_solve_err_mode[integer] :: 1: exit vel solve based on nonlin residual 2: exit based on “fixed point” metric (|u - u_last| / |u| < tol) where | | is infty-norm 3: exit based on change of norm%energysavedays[type(time_type)] :: The interval between writing the energies and other integral quantities of the run.%energysavedays_geometric[type(time_type)] :: The starting interval for computing a geometric progression of time deltas between calls to write_energy. This interval will increase by a factor of 2. after each call to write_energy.%energysave_geometric[logical] :: Logical to control whether calls to write_energy should follow a geometric progression.%write_energy_time[type(time_type)] :: The next time to write to the energy file.%geometric_end_time[type(time_type)] :: Time at which to stop the geometric progression of calls to write_energy and revert to the standard energysavedays interval.%timeunit[real] :: The length of the units for the time axis and certain input parameters including ENERGYSAVEDAYS [s].%start_time[type(time_type)] :: The start time of the simulation.%prev_is_energy_calls[integer] :: The number of times write_ice_shelf_energy has been called.%is_fileenergy_ascii[integer] :: The unit number of the ascii version of the energy file.%is_energyfile[character (len=200)] :: The name of the ice sheet energy file with path.%debug[logical] :: If true, write verbose checksums for debugging purposes and use reproducible sums.%module_is_initialized[logical] :: True if this module has been initialized.%diag[type( diag_ctrl ),pointer] :: A structure that is used to control diagnostic output.
-
type
mom_ice_shelf_dynamics/loop_bounds_type¶ A container for loop bounds.
- Type fields:
%ish[integer,private] :: Loop bounds.%ieh[integer,private] :: Loop bounds.%jsh[integer,private] :: Loop bounds.%jeh[integer,private] :: Loop bounds.
Function/Subroutine Documentation¶
-
function
mom_ice_shelf_dynamics/slope_limiter(num, denom) [real]¶ used for flux limiting in advective subroutines Van Leer limiter (source: Wikipedia) The return value is between 0 and 2 [nondim].
- Parameters:
num :: [in] The numerator of the ratio used in the Van Leer slope limiter
denom :: [in] The denominator of the ratio used in the Van Leer slope limiter
- Called from:
ice_shelf_advect_temp_xice_shelf_advect_temp_yice_shelf_advect_thickness_xice_shelf_advect_thickness_y
-
function
mom_ice_shelf_dynamics/quad_area(X, Y) [real]¶ Calculate area of quadrilateral.
- Parameters:
x :: [in] The x-positions of the vertices of the quadrilateral [L ~> m].
y :: [in] The y-positions of the vertices of the quadrilateral [L ~> m].
- Called from:
-
subroutine
mom_ice_shelf_dynamics/register_ice_shelf_dyn_restarts(G, US, param_file, CS, restart_CS)¶ This subroutine is used to register any fields related to the ice shelf dynamics that should be written to or read from the restart file.
- Parameters:
g :: [inout] The grid type describing the ice shelf grid.
us :: [in] A structure containing unit conversion factors
param_file :: [in] A structure to parse for run-time parameters
cs :: A pointer to the ice shelf dynamics control structure
restart_cs :: [inout] MOM restart control struct
- Call to:
-
subroutine
mom_ice_shelf_dynamics/initialize_ice_shelf_dyn(param_file, Time, ISS, CS, G, US, diag, new_sim, Cp_ice, Input_start_time, directory, solo_ice_sheet_in)¶ Initializes shelf model data, parameters and diagnostics.
- Parameters:
param_file :: [in] A structure to parse for run-time parameters
time :: [inout] The clock that that will indicate the model time
iss :: [in] A structure with elements that describe the ice-shelf state
cs :: A pointer to the ice shelf dynamics control structure
g :: [inout] The grid type describing the ice shelf grid.
us :: [in] A structure containing unit conversion factors
diag :: [in] A structure that is used to regulate the diagnostic output.
new_sim :: [in] If true this is a new simulation, otherwise has been started from a restart file.
cp_ice :: [in] Heat capacity of ice [Q C-1 ~> J kg-1 degC-1]
input_start_time :: [in] The start time of the simulation.
directory :: [in] The directory where the ice sheet energy file goes.
solo_ice_sheet_in :: [in] If present, this indicates whether a solo ice-sheet driver.
- Call to:
bilinear_shape_fn_gridbilinear_shape_fn_grid_1qpbilinear_shape_functions_subgridcalc_shelf_driving_stresscalc_shelf_taubcalc_shelf_viscmom_ice_shelf_initialize::initialize_ice_aglenmom_error_handler::mom_errormom_error_handler::mom_mesgupdate_od_ffrac_uncoupledupdate_velocity_masks
-
subroutine
mom_ice_shelf_dynamics/initialize_diagnostic_fields(CS, ISS, G, US, Time)¶ - Parameters:
cs :: [inout] A pointer to the ice shelf control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
time :: [in] The current model time
- Call to:
-
function
mom_ice_shelf_dynamics/ice_time_step_cfl(CS, ISS, G) [real]¶ This function returns the global maximum advective timestep that can be taken based on the current ice velocities. Because it involves finding a global minimum, it can be surprisingly expensive.
- Parameters:
cs :: [inout] The ice shelf dynamics control structure
iss :: [inout] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
- Return:
undefined :: The maximum permitted timestep based on the ice velocities [T ~> s].
-
subroutine
mom_ice_shelf_dynamics/update_ice_shelf(CS, ISS, G, US, time_step, Time, calve_ice_shelf_bergs, ocean_mass, coupled_grounding, must_update_vel)¶ This subroutine updates the ice shelf velocities, mass, stresses and properties due to the ice shelf dynamics.
- Parameters:
cs :: [inout] The ice shelf dynamics control structure
iss :: [inout] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
time_step :: [in] time step [T ~> s]
time :: [in] The current model time
calve_ice_shelf_bergs :: [in] To convert ice flux through front to bergs
ocean_mass :: [in] If present this is the mass per unit area
coupled_grounding :: [in] If true, the grounding line is determined by coupled ice-ocean dynamics
must_update_vel :: [in] Always update the ice velocities if true.
- Call to:
ice_shelf_advectice_shelf_solve_outerupdate_od_ffracupdate_od_ffrac_uncoupled
-
subroutine
mom_ice_shelf_dynamics/volume_above_floatation(CS, G, ISS, vaf, hemisphere)¶ - Parameters:
cs :: [in] The ice shelf dynamics control structure
g :: [in] The grid structure used by the ice shelf.
iss :: [in] A structure with elements that describe the ice-shelf state
vaf :: [out] area integrated volume above floatation [Z L2 ~> m3]
hemisphere :: [in] 0 for Antarctica only, 1 for Greenland only. Otherwise, all ice sheets
- Call to:
- Called from:
mom_ice_shelf::process_and_post_scalar_datamom_ice_shelf::shelf_calc_fluxmom_ice_shelf::solo_step_ice_shelf
-
subroutine
mom_ice_shelf_dynamics/masked_var_grounded(G, CS, var, varout)¶ multiplies a variable with the ice sheet grounding fraction
- Parameters:
g :: [in] The grid structure used by the ice shelf.
cs :: [in] The ice shelf dynamics control structure
var :: [in] variable in
varout :: [out] variable out
- Called from:
-
subroutine
mom_ice_shelf_dynamics/is_dynamics_post_data(time_step, Time, CS, ISS, G)¶ Ice shelf dynamics post_data calls.
- Parameters:
time_step :: Length of time for post data averaging [T ~> s].
time :: [in] The current model time
cs :: [inout] The ice shelf dynamics control structure
iss :: [inout] A structure with elements that describe the ice-shelf state
g :: [in] The grid structure used by the ice shelf.
- Call to:
mom_is_diag_mediator::disable_averagingmom_is_diag_mediator::enable_averagesice_visc_diagis_dynamics_post_data_2
-
subroutine
mom_ice_shelf_dynamics/ice_visc_diag(CS, G, ice_visc)¶ Calculate cell-centered, area-averaged, vertically integrated ice viscosity for diagnostics.
- Parameters:
cs :: [in] The ice shelf dynamics control structure
g :: [in] The grid structure used by the ice shelf.
ice_visc :: [out] area-averaged vertically integrated ice viscosity [R L2 Z T-1 ~> Pa s m]
- Called from:
-
subroutine
mom_ice_shelf_dynamics/write_ice_shelf_energy(CS, G, US, mass, area, day, time_step)¶ Writes the total ice shelf kinetic energy and mass to an ascii file.
- Parameters:
cs :: [inout] The ice shelf dynamics control structure
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
mass :: [in] The mass per unit area of the ice shelf
area :: [in] The ice shelf or ice sheet area [L2 ~> m2]
day :: [in] The current model time.
time_step :: [in] The current time step
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_advect(CS, ISS, G, time_step, Time, calve_ice_shelf_bergs)¶ This subroutine takes the velocity (on the Bgrid) and timesteps h_t = - div (uh) once. Additionally, it will update the volume of ice in partially-filled cells, and update hmask accordingly.
- Parameters:
cs :: [inout] The ice shelf dynamics control structure
iss :: [inout] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
time_step :: [in] time step [T ~> s]
time :: [in] The current model time
calve_ice_shelf_bergs :: [in] If true, track ice shelf flux through a static ice shelf, so that it can be converted into icebergs
- Call to:
calve_to_maskice_shelf_advect_thickness_xice_shelf_advect_thickness_yice_shelf_min_thickness_calvemom_error_handler::mom_errorshelf_advance_frontupdate_velocity_masks- Called from:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_solve_outer(CS, ISS, G, US, u_shlf, v_shlf, taudx, taudy, iters, Time)¶ This subroutine computes u- and v-velocities of the ice shelf iterating on non-linear ice viscosity.
- Parameters:
cs :: [inout] The ice shelf dynamics control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
u_shlf :: [inout] The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
v_shlf :: [inout] The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
iters :: [out] The number of iterations used in the solver.
time :: [in] The current model time
taudx :: [out] Driving x-stress at q-points [R L3 Z T-2 ~> kg m s-2]
taudy :: [out] Driving y-stress at q-points [R L3 Z T-2 ~> kg m s-2]
- Call to:
calc_shelf_driving_stresscalc_shelf_taubcalc_shelf_visccg_actionice_shelf_solve_innerinterpolate_h_to_bmom_error_handler::mom_mesg- Called from:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_solve_inner(CS, ISS, G, US, u_shlf, v_shlf, taudx, taudy, H_node, float_cond, hmask, conv_flag, iters, time, Phi, Phisub)¶ - Parameters:
cs :: [in] A pointer to the ice shelf control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
u_shlf :: [inout] The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
v_shlf :: [inout] The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
taudx :: [in] The x-direction driving stress [R L3 Z T-2 ~> kg m s-2]
taudy :: [in] The y-direction driving stress [R L3 Z T-2 ~> kg m s-2]
h_node :: [in] The ice shelf thickness at nodal (corner)
float_cond :: [in] If GL_regularize=true, indicates cells containing
hmask :: [in] A mask indicating which tracer points are
conv_flag :: [out] A flag indicating whether (1) or not (0) the iterations have converged to the specified tolerance
iters :: [out] The number of iterations used in the solver.
time :: [in] The current model time
phi :: [in] The gradients of bilinear basis elements at Gaussian
phisub :: [in] Quadrature structure weights at subgridscale
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_advect_thickness_x(CS, G, LB, time_step, hmask, h0, h_after_uflux, uh_ice)¶ - Parameters:
cs :: [in] A pointer to the ice shelf control structure
g :: [in] The grid structure used by the ice shelf.
lb :: [in] Loop bounds structure.
time_step :: [in] The time step for this update [T ~> s].
hmask :: [inout] A mask indicating which tracer points are
h0 :: [in] The initial ice shelf thicknesses [Z ~> m].
h_after_uflux :: [inout] The ice shelf thicknesses after
uh_ice :: [inout] The accumulated zonal ice volume flux [Z L2 ~> m3]
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_advect_thickness_y(CS, G, LB, time_step, hmask, h0, h_after_vflux, vh_ice)¶ - Parameters:
cs :: [in] A pointer to the ice shelf control structure
g :: [in] The grid structure used by the ice shelf.
lb :: [in] Loop bounds structure.
time_step :: [in] The time step for this update [T ~> s].
hmask :: [inout] A mask indicating which tracer points are
h0 :: [in] The initial ice shelf thicknesses [Z ~> m].
h_after_vflux :: [inout] The ice shelf thicknesses after
vh_ice :: [inout] The accumulated meridional ice volume flux [Z L2 ~> m3]
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/shelf_advance_front(CS, ISS, G, hmask, uh_ice, vh_ice)¶ - Parameters:
cs :: [in] A pointer to the ice shelf control structure
iss :: [inout] A structure with elements that describe the ice-shelf state
g :: [in] The grid structure used by the ice shelf.
hmask :: [inout] A mask indicating which tracer points are
uh_ice :: [inout] The accumulated zonal ice volume flux [Z L2 ~> m3]
vh_ice :: [inout] The accumulated meridional ice volume flux [Z L2 ~> m3]
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_min_thickness_calve(G, h_shelf, area_shelf_h, hmask, thickness_calve, halo)¶ Apply a very simple calving law using a minimum thickness rule.
- Parameters:
g :: [in] The grid structure used by the ice shelf.
h_shelf :: [inout] The ice shelf thickness [Z ~> m].
area_shelf_h :: [inout] The area per cell covered by the ice shelf [L2 ~> m2].
hmask :: [inout] A mask indicating which tracer points are partly or fully covered by an ice-shelf
thickness_calve :: [in] The thickness at which to trigger calving [Z ~> m].
halo :: [in] The number of halo points to use. If not present, work on the entire data domain.
- Called from:
-
subroutine
mom_ice_shelf_dynamics/calve_to_mask(G, h_shelf, area_shelf_h, hmask, calve_mask)¶ - Parameters:
g :: [in] The grid structure used by the ice shelf.
h_shelf :: [inout] The ice shelf thickness [Z ~> m].
area_shelf_h :: [inout] The area per cell covered by the ice shelf [L2 ~> m2].
hmask :: [inout] A mask indicating which tracer points are partly or fully covered by an ice-shelf
calve_mask :: [in] A mask that indicates where the ice shelf can exist, and where it will calve.
- Called from:
-
subroutine
mom_ice_shelf_dynamics/calc_shelf_driving_stress(CS, ISS, G, US, taudx, taudy, OD)¶ Calculate driving stress using cell-centered bed elevation and ice thickness.
- Parameters:
cs :: [in] A pointer to the ice shelf control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
od :: [in] ocean floor depth at tracer points [Z ~> m].
taudx :: [inout] X-direction driving stress at q-points [R L3 Z T-2 ~> kg m s-2]
taudy :: [inout] Y-direction driving stress at q-points [R L3 Z T-2 ~> kg m s-2]
- Called from:
-
subroutine
mom_ice_shelf_dynamics/cg_action(CS, uret, vret, u_shlf, v_shlf, Phi, Phisub, umask, vmask, hmask, H_node, ice_visc, float_cond, bathyT, basal_trac, G, US, is, ie, js, je, dens_ratio)¶ - Parameters:
cs :: [in] A pointer to the ice shelf control structure
g :: [in] The grid structure used by the ice shelf.
uret :: [inout] The retarding stresses working at u-points [R L3 Z T-2 ~> kg m s-2].
vret :: [inout] The retarding stresses working at v-points [R L3 Z T-2 ~> kg m s-2].
phi :: [in] The gradients of bilinear basis elements at Gaussian
phisub :: [in] Quadrature structure weights at subgridscale
u_shlf :: [in] The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
v_shlf :: [in] The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
umask :: [in] A coded mask indicating the nature of the
vmask :: [in] A coded mask indicating the nature of the
h_node :: [in] The ice shelf thickness at nodal (corner)
hmask :: [in] A mask indicating which tracer points are
ice_visc :: [in] A field related to the ice viscosity from Glen’s
float_cond :: [in] If GL_regularize=true, an array indicating where the ice
bathyt :: [in] The depth of ocean bathymetry at tracer points
basal_trac :: [in] Area-integrated taub_beta field related to the nonlinear
dens_ratio :: [in] The density of ice divided by the density of seawater, nondimensional
us :: [in] A structure containing unit conversion factors
is :: [in] The starting i-index to work on
ie :: [in] The ending i-index to work on
js :: [in] The starting j-index to work on
je :: [in] The ending j-index to work on
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/cg_action_subgrid_basal(Phisub, H, U, V, bathyT, dens_ratio, Ucontr, Vcontr)¶ - Parameters:
phisub :: [in] Quadrature structure weights at subgridscale
h :: [in] The ice shelf thickness at nodal (corner) points [Z ~> m].
u :: [in] The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
v :: [in] The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
bathyt :: [in] The depth of ocean bathymetry at tracer points relative to sea-level [Z ~> m].
dens_ratio :: [in] The density of ice divided by the density of seawater [nondim]
ucontr :: [out] The areal average of u-velocities where the ice shelf is grounded, or 0 where it is floating [L T-1 ~> m s-1].
vcontr :: [out] The areal average of v-velocities where the ice shelf is grounded, or 0 where it is floating [L T-1 ~> m s-1].
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/sum_square_matrix(sum_out, mat_in, n)¶ - Parameters:
n :: [in] The length and width of each matrix in mat_in
mat_in :: [in] The n x n matrix whose elements will be summed
sum_out :: [out] The sum of the elements of matrix mat_in
- Called from:
-
subroutine
mom_ice_shelf_dynamics/matrix_diagonal(CS, G, US, float_cond, H_node, ice_visc, basal_trac, hmask, dens_ratio, Phi, Phisub, u_diagonal, v_diagonal)¶ returns the diagonal entries of the matrix for a Jacobi preconditioning
- Parameters:
cs :: [in] A pointer to the ice shelf control structure
g :: [in] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
float_cond :: [in] If GL_regularize=true, indicates cells containing
h_node :: [in] The ice shelf thickness at nodal
ice_visc :: [in] A field related to the ice viscosity from Glen’s
basal_trac :: [in] Area-integrated taub_beta field related to the nonlinear
hmask :: [in] A mask indicating which tracer points are
dens_ratio :: [in] The density of ice divided by the density of seawater [nondim]
phi :: [in] The gradients of bilinear basis elements at Gaussian
phisub :: [in] Quadrature structure weights at subgridscale locations for finite element calculations [nondim]
u_diagonal :: [inout] The diagonal elements of the u-velocity
v_diagonal :: [inout] The diagonal elements of the v-velocity
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/cg_diagonal_subgrid_basal(Phisub, H_node, bathyT, dens_ratio, f_grnd)¶ - Parameters:
phisub :: [in] Quadrature structure weights at subgridscale
h_node :: [in] The ice shelf thickness at nodal (corner) points [Z ~> m].
bathyt :: [in] The depth of ocean bathymetry at tracer points [Z ~> m].
dens_ratio :: [in] The density of ice divided by the density of seawater [nondim]
f_grnd :: [out] The weighted fraction of the sub-cell where the ice shelf is grounded [nondim]
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/is_dynamics_post_data_2(CS, ISS, G)¶ Post_data calls related to ice-sheet flux divergence, strain-rate, and deviatoric stress.
- Parameters:
cs :: [inout] A pointer to the ice shelf control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [in] The grid structure used by the ice shelf.
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/calc_shelf_visc(CS, ISS, G, US, u_shlf, v_shlf)¶ Update depth integrated viscosity, based on horizontal strain rates.
- Parameters:
cs :: [inout] A pointer to the ice shelf control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [in] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
u_shlf :: [inout] The zonal ice shelf velocity [L T-1 ~> m s-1].
v_shlf :: [inout] The meridional ice shelf velocity [L T-1 ~> m s-1].
- Called from:
-
subroutine
mom_ice_shelf_dynamics/calc_shelf_taub(CS, ISS, G, US, u_shlf, v_shlf)¶ Update basal shear.
- Parameters:
cs :: [inout] A pointer to the ice shelf control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [in] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
u_shlf :: [inout] The zonal ice shelf velocity [L T-1 ~> m s-1].
v_shlf :: [inout] The meridional ice shelf velocity [L T-1 ~> m s-1].
- Called from:
-
subroutine
mom_ice_shelf_dynamics/update_od_ffrac(CS, G, US, ocean_mass, find_avg)¶ - Parameters:
cs :: [inout] A pointer to the ice shelf control structure
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
ocean_mass :: [in] The mass per unit area of the ocean [R Z ~> kg m-2].
find_avg :: [in] If true, find the average of OD and ffrac, and reset the underlying running sums to 0.
- Called from:
-
subroutine
mom_ice_shelf_dynamics/update_od_ffrac_uncoupled(CS, G, h_shelf)¶ - Parameters:
cs :: [inout] A pointer to the ice shelf control structure
g :: [in] The grid structure used by the ice shelf.
h_shelf :: [in] the thickness of the ice shelf [Z ~> m].
- Called from:
-
subroutine
mom_ice_shelf_dynamics/change_in_draft(CS, G, h_shelf0, h_shelf1, ddraft)¶ - Parameters:
cs :: [inout] A pointer to the ice shelf control structure
g :: [in] The grid structure used by the ice shelf.
h_shelf0 :: [in] the previous thickness of the ice shelf [Z ~> m].
h_shelf1 :: [in] the current thickness of the ice shelf [Z ~> m].
ddraft :: [inout] the change in shelf draft thickness
-
subroutine
mom_ice_shelf_dynamics/bilinear_shape_functions(X, Y, Phi, area)¶ This subroutine calculates the gradients of bilinear basis elements that that are centered at the vertices of the cell. Values are calculated at points of gaussian quadrature.
- Parameters:
x :: [in] The x-positions of the vertices of the quadrilateral [L ~> m].
y :: [in] The y-positions of the vertices of the quadrilateral [L ~> m].
phi :: [inout] The gradients of bilinear basis elements at Gaussian quadrature points surrounding the cell vertices [L-1 ~> m-1].
area :: [out] The quadrilateral cell area [L2 ~> m2].
- Call to:
-
subroutine
mom_ice_shelf_dynamics/bilinear_shape_fn_grid(G, i, j, Phi)¶ This subroutine calculates the gradients of bilinear basis elements that are centered at the vertices of the cell using a locally orthogoal MOM6 grid. Values are calculated at points of gaussian quadrature.
- Parameters:
g :: [in] The grid structure used by the ice shelf.
i :: [in] The i-index in the grid to work on.
j :: [in] The j-index in the grid to work on.
phi :: [inout] The gradients of bilinear basis elements at Gaussian quadrature points surrounding the cell vertices [L-1 ~> m-1].
- Called from:
-
subroutine
mom_ice_shelf_dynamics/bilinear_shape_fn_grid_1qp(G, i, j, Phi)¶ This subroutine calculates the gradients of bilinear basis elements that are centered at the vertices of the cell using a locally orthogoal MOM6 grid. Values are calculated at a sinlge cell-centered quadrature point, which should match the grid cell h-point.
- Parameters:
g :: [in] The grid structure used by the ice shelf.
i :: [in] The i-index in the grid to work on.
j :: [in] The j-index in the grid to work on.
phi :: [inout] The gradients of bilinear basis elements at Gaussian quadrature points surrounding the cell vertices [L-1 ~> m-1].
- Called from:
-
subroutine
mom_ice_shelf_dynamics/bilinear_shape_functions_subgrid(Phisub, nsub)¶ - Parameters:
nsub :: [in] The number of subgridscale quadrature locations in each direction
phisub :: [inout] Quadrature structure weights at subgridscale
- Called from:
-
subroutine
mom_ice_shelf_dynamics/update_velocity_masks(CS, G, hmask, umask, vmask, u_face_mask, v_face_mask)¶ - Parameters:
cs :: [in] A pointer to the ice shelf dynamics control structure
g :: [inout] The grid structure used by the ice shelf.
hmask :: [in] A mask indicating which tracer points are
umask :: [out] A coded mask indicating the nature of the
vmask :: [out] A coded mask indicating the nature of the
u_face_mask :: [out] A coded mask for velocities at the C-grid u-face
v_face_mask :: [out] A coded mask for velocities at the C-grid v-face
- Called from:
-
subroutine
mom_ice_shelf_dynamics/interpolate_h_to_b(G, h_shelf, hmask, H_node, min_h_shelf)¶ Interpolate the ice shelf thickness from tracer point to nodal points, subject to a mask.
- Parameters:
g :: [in] The grid structure used by the ice shelf.
h_shelf :: [in] The ice shelf thickness at tracer points [Z ~> m].
hmask :: [in] A mask indicating which tracer points are
h_node :: [inout] The ice shelf thickness at nodal (corner)
min_h_shelf :: [in] The minimum ice thickness used during ice dynamics [Z ~> m].
- Called from:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_dyn_end(CS)¶ Deallocates all memory associated with the ice shelf dynamics module.
- Parameters:
cs :: A pointer to the ice shelf dynamics control structure
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_temp(CS, ISS, G, US, time_step, melt_rate, Time)¶ This subroutine updates the vertically averaged ice shelf temperature.
- Parameters:
cs :: [inout] A pointer to the ice shelf control structure
iss :: [in] A structure with elements that describe the ice-shelf state
g :: [inout] The grid structure used by the ice shelf.
us :: [in] A structure containing unit conversion factors
time_step :: [in] The time step for this update [T ~> s].
melt_rate :: [in] basal melt rate [R Z T-1 ~> kg m-2 s-1]
time :: [in] The current model time
- Call to:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_advect_temp_x(CS, G, time_step, hmask, h0, h_after_uflux)¶ - Parameters:
cs :: [in] A pointer to the ice shelf control structure
g :: [inout] The grid structure used by the ice shelf.
time_step :: [in] The time step for this update [T ~> s].
hmask :: [in] A mask indicating which tracer points are
h0 :: [in] The initial ice shelf thicknesses times temperature [C Z ~> degC m]
h_after_uflux :: [inout] The ice shelf thicknesses times temperature after
- Call to:
- Called from:
-
subroutine
mom_ice_shelf_dynamics/ice_shelf_advect_temp_y(CS, G, time_step, hmask, h_after_uflux, h_after_vflux)¶ - Parameters:
cs :: [in] A pointer to the ice shelf control structure
g :: [in] The grid structure used by the ice shelf.
time_step :: [in] The time step for this update [T ~> s].
hmask :: [in] A mask indicating which tracer points are
h_after_uflux :: [in] The ice shelf thicknesses times temperature after
h_after_vflux :: [inout] The ice shelf thicknesses times temperature after
- Call to:
- Called from: