MOM_ice_shelf_dynamics.F90
1! This file is part of MOM6, the Modular Ocean Model version 6.
2! See the LICENSE file for licensing information.
3! SPDX-License-Identifier: Apache-2.0
4
5!> Implements a crude placeholder for a later implementation of full
6!! ice shelf dynamics.
7module mom_ice_shelf_dynamics
8
9use mom_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
10use mom_cpu_clock, only : clock_component, clock_routine
11use mom_is_diag_mediator, only : post_data=>post_is_data
12use mom_is_diag_mediator, only : register_diag_field=>register_mom_is_diag_field, safe_alloc_ptr
13!use MOM_IS_diag_mediator, only : MOM_IS_diag_mediator_init, set_IS_diag_mediator_grid
14use mom_is_diag_mediator, only : diag_ctrl, time_type, enable_averages, disable_averaging
15use mom_domains, only : mom_domains_init, clone_mom_domain
16use mom_domains, only : pass_var, pass_vector, to_all, cgrid_ne, bgrid_ne, agrid, corner, center
17use mom_error_handler, only : mom_error, mom_mesg, fatal, warning, is_root_pe
18use mom_file_parser, only : read_param, get_param, log_param, log_version, param_file_type
19use mom_grid, only : mom_grid_init, ocean_grid_type
20use mom_io, only : file_exists, slasher, mom_read_data
21use mom_io, only : open_ascii_file, get_filename_appendix
22use mom_io, only : append_file, writeonly_file
23use mom_restart, only : register_restart_field, mom_restart_cs
24use mom_time_manager, only : time_type, get_time, set_time, time_type_to_real, operator(>)
25use mom_time_manager, only : operator(+), operator(-), operator(*), operator(/)
26use mom_time_manager, only : operator(/=), operator(<=), operator(>=), operator(<)
27use mom_unit_scaling, only : unit_scale_type, unit_scaling_init
28!MJH use MOM_ice_shelf_initialize, only : initialize_ice_shelf_boundary
29use mom_ice_shelf_state, only : ice_shelf_state
30use mom_coms, only : reproducing_sum, max_across_pes, min_across_pes
31use mom_checksums, only : hchksum, qchksum
32use mom_ice_shelf_initialize, only : initialize_ice_shelf_boundary_channel,initialize_ice_flow_from_file
33use mom_ice_shelf_initialize, only : initialize_ice_shelf_boundary_from_file,initialize_ice_c_basal_friction
34use mom_ice_shelf_initialize, only : initialize_ice_aglen
35implicit none ; private
36
37#include <MOM_memory.h>
38
39public register_ice_shelf_dyn_restarts, initialize_ice_shelf_dyn, update_ice_shelf, is_dynamics_post_data
42public masked_var_grounded
43
44! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
45! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
46! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
47! vary with the Boussinesq approximation, the Boussinesq variant is given first.
48
49!> The control structure for the ice shelf dynamics.
50type, public :: ice_shelf_dyn_cs ; private
51 real, pointer, dimension(:,:) :: u_shelf => null() !< the zonal velocity of the ice shelf/sheet
52 !! on q-points (B grid) [L T-1 ~> m s-1]
53 real, pointer, dimension(:,:) :: v_shelf => null() !< the meridional velocity of the ice shelf/sheet
54 !! on q-points (B grid) [L T-1 ~> m s-1]
55 real, pointer, dimension(:,:) :: taudx_shelf => null() !< the zonal driving stress of the ice shelf/sheet
56 !! on q-points (C grid) [R L2 T-2 ~> Pa]
57 real, pointer, dimension(:,:) :: taudy_shelf => null() !< the meridional driving stress of the ice shelf/sheet
58 !! on q-points (C grid) [R L2 T-2 ~> Pa]
59 real, pointer, dimension(:,:) :: sx_shelf => null() !< the zonal surface slope of the ice shelf/sheet
60 !! on q-points (B grid) [nondim]
61 real, pointer, dimension(:,:) :: sy_shelf => null() !< the meridional surface slope of the ice shelf/sheet
62 !! on q-points (B grid) [nondim]
63 real, pointer, dimension(:,:) :: u_face_mask => null() !< mask for velocity boundary conditions on the C-grid
64 !! u-face - this is because the FEM cares about FACES THAT GET INTEGRATED OVER,
65 !! not vertices. Will represent boundary conditions on computational boundary
66 !! (or permanent boundary between fast-moving and near-stagnant ice
67 !! FOR NOW: 1=interior bdry, 0=no-flow boundary, 2=stress bdry condition,
68 !! 3=inhomogeneous Dirichlet boundary for u and v, 4=flux boundary: at these
69 !! faces a flux will be specified which will override velocities; a homogeneous
70 !! velocity condition will be specified (this seems to give the solver less
71 !! difficulty) 5=inhomogenous Dirichlet boundary for u only. 6=inhomogenous
72 !! Dirichlet boundary for v only
73 real, pointer, dimension(:,:) :: v_face_mask => null() !< A mask for velocity boundary conditions on the C-grid
74 !! v-face, with valued defined similarly to u_face_mask, but 5 is Dirichlet for v
75 !! and 6 is Dirichlet for u
76 real, pointer, dimension(:,:) :: u_face_mask_bdry => null() !< A duplicate copy of u_face_mask?
77 real, pointer, dimension(:,:) :: v_face_mask_bdry => null() !< A duplicate copy of v_face_mask?
78 real, pointer, dimension(:,:) :: u_flux_bdry_val => null() !< The ice volume flux per unit face length into the cell
79 !! through open boundary u-faces (where u_face_mask=4) [Z L T-1 ~> m2 s-1]
80 real, pointer, dimension(:,:) :: v_flux_bdry_val => null() !< The ice volume flux per unit face length into the cell
81 !! through open boundary v-faces (where v_face_mask=4) [Z L T-1 ~> m2 s-1]??
82 ! needed where u_face_mask is equal to 4, similarly for v_face_mask
83 real, pointer, dimension(:,:) :: umask => null() !< u-mask on the actual degrees of freedom (B grid)
84 !! 1=normal node, 3=inhomogeneous boundary node,
85 !! 0 - no flow node (will also get ice-free nodes)
86 real, pointer, dimension(:,:) :: vmask => null() !< v-mask on the actual degrees of freedom (B grid)
87 !! 1=normal node, 3=inhomogeneous boundary node,
88 !! 0 - no flow node (will also get ice-free nodes)
89 real, pointer, dimension(:,:) :: calve_mask => null() !< a mask to prevent the ice shelf front from
90 !! advancing past its initial position (but it may retreat)
91 real, pointer, dimension(:,:) :: t_shelf => null() !< Vertically integrated temperature in the ice shelf/stream,
92 !! on corner-points (B grid) [C ~> degC]
93 real, pointer, dimension(:,:) :: tmask => null() !< A mask on tracer points that is 1 where there is ice.
94 real, pointer, dimension(:,:,:) :: ice_visc => null() !< Area and depth-integrated Glen's law ice viscosity
95 !! (Pa m3 s) in [R L4 Z T-1 ~> kg m2 s-1].
96 !! at either 1 (cell-centered) or 4 quadrature points per cell
97 real, pointer, dimension(:,:,:) :: newton_visc_factor => null() !< Newton tangent stiffness coefficient:
98 !! (1/n_glen - 1)/2 * ice_visc / eps_e2 at each
99 !! viscosity quadrature point [R L4 Z T ~> kg m2 s]
100 real, pointer, dimension(:,:,:) :: newton_str_ux => null() !< Longitudinal x-strain-rate ux at each viscosity
101 !! quadrature point for Newton iterations [T-1 ~> s-1]
102 real, pointer, dimension(:,:,:) :: newton_str_vy => null() !< Longitudinal y-strain-rate vy at each viscosity
103 !! quadrature point for Newton iterations [T-1 ~> s-1]
104 real, pointer, dimension(:,:,:) :: newton_str_sh => null() !< Engineering shear strain-rate uy+vx at each
105 !! viscosity quadrature point for Newton iterations [T-1 ~> s-1]
106 real, pointer, dimension(:,:) :: newton_umid => null() !< Cell-averaged zonal velocity u at the current outer
107 !! iterate, for Newton basal drag correction [L T-1 ~> m s-1]
108 real, pointer, dimension(:,:) :: newton_vmid => null() !< Cell-averaged meridional velocity v at the current
109 !! outer iterate, for Newton basal drag correction [L T-1 ~> m s-1]
110 real, pointer, dimension(:,:) :: newton_drag_coef => null() !< Newton basal drag correction coefficient:
111 !! 2 * d(basal_trac)/d(|u|^2) * area = d(tau_b_i)/d(u_j) - basal_trac*delta_ij
112 !! expressed as the u_i*u_j tensor coefficient [R Z T ~> kg m-2 s]
113 real, pointer, dimension(:,:) :: aglen_visc => null() !< Ice-stiffness parameter in Glen's law ice viscosity,
114 !! often in [Pa-3 s-1] if n_Glen is 3.
115 real, pointer, dimension(:,:) :: u_bdry_val => null() !< The zonal ice velocity at inflowing boundaries
116 !! [L yr-1 ~> m yr-1]
117 real, pointer, dimension(:,:) :: v_bdry_val => null() !< The meridional ice velocity at inflowing boundaries
118 !! [L yr-1 ~> m yr-1]
119 real, pointer, dimension(:,:) :: h_bdry_val => null() !< The ice thickness at inflowing boundaries [Z ~> m].
120 real, pointer, dimension(:,:) :: t_bdry_val => null() !< The ice temperature at inflowing boundaries [C ~> degC].
121
122 real, pointer, dimension(:,:) :: bed_elev => null() !< The bed elevation used for ice dynamics [Z ~> m],
123 !! relative to mean sea-level. This is
124 !! the same as G%bathyT+Z_ref, when below sea-level.
125 !! Sign convention: positive below sea-level, negative above.
126
127 real, pointer, dimension(:,:) :: basal_traction => null() !< The area-integrated taub_beta field
128 !! (m2 Pa s m-1, or kg s-1) related to the nonlinear part
129 !! of "linearized" basal stress (Pa) [R Z L2 T-1 ~> kg s-1]
130 !! The exact form depends on basal law exponent and/or whether flow is "hybridized" a la Goldberg 2011
131 real, pointer, dimension(:,:) :: c_basal_friction => null()!< Coefficient in sliding law tau_b = C u^(n_basal_fric),
132 !! units of [R L Z T-2 (s m-1)^(n_basal_fric) ~> Pa (s m-1)^(n_basal_fric)]
133 real, pointer, dimension(:,:) :: od_rt => null() !< A running total for calculating OD_av [Z ~> m].
134 real, pointer, dimension(:,:) :: ground_frac_rt => null() !< A running total for calculating ground_frac.
135 real, pointer, dimension(:,:) :: od_av => null() !< The time average open ocean depth [Z ~> m].
136 real, pointer, dimension(:,:) :: ground_frac => null() !< Fraction of the time a cell is "exposed", i.e. the column
137 !! thickness is below a threshold and interacting with the rock [nondim]. When this
138 !! is 1, the ice-shelf is grounded
139 real, pointer, dimension(:,:) :: float_cond => null() !< If GL_regularize=true, indicates cells containing
140 !! the grounding line (float_cond=1) or not (float_cond=0)
141 real, pointer, dimension(:,:,:,:) :: phi => null() !< The gradients of bilinear basis elements at Gaussian
142 !! 4 quadrature points surrounding the cell vertices [L-1 ~> m-1].
143 real, pointer, dimension(:,:,:) :: phic => null() !< The gradients of bilinear basis elements at 1 cell-centered
144 !! quadrature point per cell [L-1 ~> m-1].
145 real, pointer, dimension(:,:,:,:,:,:) :: phisub => null() !< Quadrature structure weights at subgridscale
146 !! locations for finite element calculations [nondim]
147 integer :: od_rt_counter = 0 !< A counter of the number of contributions to OD_rt.
148
149 real :: velocity_update_time_step !< The time interval over which to update the ice shelf velocity
150 !! using the nonlinear elliptic equation, or 0 to update every timestep [T ~> s].
151 ! DNGoldberg thinks this should be done no more often than about once a day
152 ! (maybe longer) because it will depend on ocean values that are averaged over
153 ! this time interval, and solving for the equilibrated flow will begin to lose
154 ! meaning if it is done too frequently.
155 real :: elapsed_velocity_time !< The elapsed time since the ice velocities were last updated [T ~> s].
156
157 real :: g_earth !< The gravitational acceleration [L2 Z-1 T-2 ~> m s-2].
158 real :: density_ice !< A typical density of ice [R ~> kg m-3].
159 real :: cp_ice !< The heat capacity of fresh ice [Q C-1 ~> J kg-1 degC-1].
160
161 logical :: advect_shelf !< If true (default), advect ice shelf and evolve thickness
162 logical :: reentrant_x !< If true, the domain is zonally reentrant
163 logical :: reentrant_y !< If true, the domain is meridionally reentrant
164 logical :: alternate_first_direction_is !< If true, alternate whether the x- or y-direction
165 !! updates occur first in directionally split parts of the calculation.
166 integer :: first_direction_is !< An integer that indicates which direction is
167 !! to be updated first in directionally split
168 !! parts of the ice sheet calculation (e.g. advection).
169 real :: first_dir_restart_is = -1.0 !< A real copy of CS%first_direction_IS for use in restart files
170 integer :: visc_qps !< The number of quadrature points per cell (1 or 4) on which to calculate ice viscosity.
171 character(len=40) :: ice_viscosity_compute !< Specifies whether the ice viscosity is computed internally
172 !! according to Glen's flow law; is constant (for debugging purposes)
173 !! or using observed strain rates and read from a file
174 logical :: shelf_top_slope_bugs !< If true, use directionally inconsistent estimates of the grid
175 !! spacing when calculating the ice shelf surface slope, and underestimate
176 !! slopes near the edge of the ice shelf by a factor of 2.
177 logical :: gl_regularize !< Specifies whether to regularize the floatation condition
178 !! at the grounding line as in Goldberg Holland Schoof 2009
179 integer :: n_sub_regularize
180 !< partition of cell over which to integrate for
181 !! interpolated grounding line the (rectangular) is
182 !! divided into nxn equally-sized rectangles, over which
183 !! basal contribution is integrated (iterative quadrature)
184 logical :: gl_couple !< whether to let the floatation condition be
185 !! determined by ocean column thickness means update_OD_ffrac
186 !! will be called (note: GL_regularize and GL_couple
187 !! should be exclusive)
188
189 real :: cfl_factor !< A factor used to limit subcycled advective timestep in uncoupled runs
190 !! i.e. dt <= CFL_factor * min(dx / u) [nondim]
191
192 real :: min_h_shelf !< The minimum ice thickness used during ice dynamics [Z ~> m].
193 real :: min_basal_traction !< The minimum basal traction for grounded ice (Pa m-1 s) [R Z T-1 ~> kg m-2 s-1]
194 real :: max_surface_slope !< The maximum allowed ice-sheet surface slope (to ignore, set to zero) [nondim]
195 real :: min_ice_visc !< The minimum allowed Glen's law ice viscosity (Pa s), in [R L2 T-1 ~> kg m-1 s-1].
196
197 real :: n_glen !< Nonlinearity exponent in Glen's Law [nondim]
198 real :: eps_glen_min !< Min. strain rate to avoid infinite Glen's law viscosity, [T-1 ~> s-1].
199 real :: n_basal_fric !< Exponent in sliding law tau_b = C u^(m_slide) [nondim]
200 logical :: coulombfriction !< Use Coulomb friction law (Schoof 2005, Gagliardini et al 2007)
201 real :: cf_minn !< Minimum Coulomb friction effective pressure [R Z L T-2 ~> Pa]
202 real :: cf_postpeak !< Coulomb friction post peak exponent [nondim]
203 real :: cf_max !< Coulomb friction maximum coefficient [nondim]
204 real :: density_ocean_avg !< A typical ocean density [R ~> kg m-3]. This does not affect ocean
205 !! circulation or thermodynamics. It is used to estimate the
206 !! gravitational driving force at the shelf front (until we think of
207 !! a better way to do it, but any difference will be negligible).
208 real :: thresh_float_col_depth !< The water column depth over which the shelf if considered to be floating
209 logical :: moving_shelf_front !< Specify whether to advance shelf front (and calve).
210 logical :: calve_to_mask !< If true, calve off the ice shelf when it passes the edge of a mask.
211 real :: min_thickness_simple_calve !< min. ice shelf thickness criteria for calving [Z ~> m].
212 real :: t_shelf_missing !< An ice shelf temperature to use where there is no ice shelf [C ~> degC]
213 real :: cg_tolerance !< The tolerance in the CG solver, relative to initial residual, that
214 !! determines when to stop the conjugate gradient iterations [nondim].
215 real :: cg_tol_newton !< Working CG tolerance for the current inner solve [nondim].
216 real :: nonlinear_tolerance !< The fractional nonlinear tolerance, relative to the initial error,
217 !! that sets when to stop the iterative velocity solver [nondim]
218 real :: newton_after_tolerance !< The fractional nonlinear tolerance, relative to the initial error, at
219 !! which to switch from Picard to Newton iterations in the velocity solver [nondim]
220 logical :: newton_adapt_cg_tol !< Use an adaptive CG tolerance during Newton iterations
221 integer :: cg_max_iterations !< The maximum number of iterations that can be used in the CG solver
222 integer :: nonlin_solve_err_mode !< 1: exit vel solve based on nonlin residual
223 !! 2: exit based on "fixed point" metric (|u - u_last| / |u| < tol) where | | is infty-norm
224 !! 3: exit based on change of norm
225
226 ! for write_ice_shelf_energy
227 type(time_type) :: energysavedays !< The interval between writing the energies
228 !! and other integral quantities of the run.
229 type(time_type) :: energysavedays_geometric !< The starting interval for computing a geometric
230 !! progression of time deltas between calls to
231 !! write_energy. This interval will increase by a factor of 2.
232 !! after each call to write_energy.
233 logical :: energysave_geometric !< Logical to control whether calls to write_energy should
234 !! follow a geometric progression
235 type(time_type) :: write_energy_time !< The next time to write to the energy file.
236 type(time_type) :: geometric_end_time !< Time at which to stop the geometric progression
237 !! of calls to write_energy and revert to the standard
238 !! energysavedays interval
239 real :: timeunit !< The length of the units for the time axis and certain input parameters
240 !! including ENERGYSAVEDAYS [s].
241 type(time_type) :: start_time !< The start time of the simulation.
242 ! Start_time is set in MOM_initialization.F90
243 integer :: prev_is_energy_calls = 0 !< The number of times write_ice_shelf_energy has been called.
244 integer :: is_fileenergy_ascii !< The unit number of the ascii version of the energy file.
245 character(len=200) :: is_energyfile !< The name of the ice sheet energy file with path.
246
247 ! ids for outputting intermediate thickness in advection subroutine (debugging)
248 !integer :: id_h_after_uflux = -1, id_h_after_vflux = -1, id_h_after_adv = -1
249
250 logical :: debug !< If true, write verbose checksums for debugging purposes
251 !! and use reproducible sums
252 logical :: doing_newton = .false. !< If true, the outer iteration is using Newton (tangent) linearization
253 !! instead of Picard (secant) linearization for the ice viscosity
254 logical :: module_is_initialized = .false. !< True if this module has been initialized.
255
256 !>@{ Diagnostic handles
257 integer :: id_u_shelf = -1, id_v_shelf = -1, id_shelf_speed, id_t_shelf = -1, &
258 id_taudx_shelf = -1, id_taudy_shelf = -1, id_taud_shelf = -1, id_bed_elev = -1, &
259 id_ground_frac = -1, id_col_thick = -1, id_od_av = -1, id_float_cond = -1, &
260 id_u_mask = -1, id_v_mask = -1, id_ufb_mask =-1, id_vfb_mask = -1, id_t_mask = -1, &
261 id_sx_shelf = -1, id_sy_shelf = -1, id_surf_slope_mag_shelf, &
262 id_duhdx = -1, id_dvhdy = -1, id_fluxdiv = -1, &
263 id_strainrate_xx = -1, id_strainrate_yy = -1, id_strainrate_xy = -1, &
264 id_pstrainrate_1 = -1, id_pstrainrate_2, &
265 id_devstress_xx = -1, id_devstress_yy = -1, id_devstress_xy = -1, &
266 id_pdevstress_1 = -1, id_pdevstress_2 = -1
267
268 !>@}
269 ! ids for outputting intermediate thickness in advection subroutine (debugging)
270 !>@{ Diagnostic handles for debugging
271 integer :: id_h_after_uflux = -1, id_h_after_vflux = -1, id_h_after_adv = -1, &
272 id_visc_shelf = -1, id_taub = -1
273 !>@}
274 type(diag_ctrl), pointer :: diag => null() !< A structure that is used to control diagnostic output.
275
276end type ice_shelf_dyn_cs
277
278!> A container for loop bounds
279type :: loop_bounds_type ; private
280 integer :: ish !< Starting i-index of the computational domain [nondim]
281 integer :: ieh !< Ending i-index of the computational domain [nondim]
282 integer :: jsh !< Starting j-index of the computational domain [nondim]
283 integer :: jeh !< Ending j-index of the computational domain [nondim]
284end type loop_bounds_type
285
286contains
287
288!> used for flux limiting in advective subroutines Van Leer limiter (source: Wikipedia)
289!! The return value is between 0 and 2 [nondim].
290function slope_limiter(num, denom)
291 real, intent(in) :: num !< The numerator of the ratio used in the Van Leer slope limiter
292 real, intent(in) :: denom !< The denominator of the ratio used in the Van Leer slope limiter
293 real :: slope_limiter ! The slope limiter value, between 0 and 2 [nondim].
294 real :: r ! The ratio of num/denom [nondim]
295
296 if (denom == 0) then
297 slope_limiter = 0
298 elseif (num*denom <= 0) then
299 slope_limiter = 0
300 else
301 r = num/denom
302 slope_limiter = (r+abs(r))/(1+abs(r))
303 endif
304
305end function slope_limiter
306
307!> Calculate area of quadrilateral.
308function quad_area (X, Y)
309 real, dimension(4), intent(in) :: x !< The x-positions of the vertices of the quadrilateral [L ~> m].
310 real, dimension(4), intent(in) :: y !< The y-positions of the vertices of the quadrilateral [L ~> m].
311 real :: quad_area ! Computed area [L2 ~> m2]
312 real :: p2, q2, a2, c2, b2, d2
313
314! X and Y must be passed in the form
315 ! 3 - 4
316 ! | |
317 ! 1 - 2
318
319 p2 = ( ((x(4)-x(1))**2) + ((y(4)-y(1))**2) ) ; q2 = ( ((x(3)-x(2))**2) + ((y(3)-y(2))**2) )
320 a2 = ( ((x(3)-x(4))**2) + ((y(3)-y(4))**2) ) ; c2 = ( ((x(1)-x(2))**2) + ((y(1)-y(2))**2) )
321 b2 = ( ((x(2)-x(4))**2) + ((y(2)-y(4))**2) ) ; d2 = ( ((x(3)-x(1))**2) + ((y(3)-y(1))**2) )
322 quad_area = .25 * sqrt(4*p2*q2-(b2+d2-a2-c2)**2)
323
324end function quad_area
325
326!> This subroutine is used to register any fields related to the ice shelf
327!! dynamics that should be written to or read from the restart file.
328subroutine register_ice_shelf_dyn_restarts(G, US, param_file, CS, restart_CS)
329 type(ocean_grid_type), intent(inout) :: g !< The grid type describing the ice shelf grid.
330 type(unit_scale_type), intent(in) :: us !< A structure containing unit conversion factors
331 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
332 type(ice_shelf_dyn_cs), pointer :: cs !< A pointer to the ice shelf dynamics control structure
333 type(mom_restart_cs), intent(inout) :: restart_cs !< MOM restart control struct
334
335 ! Local variables
336 real :: t_shelf_missing ! An ice shelf temperature to use where there is no ice shelf [C ~> degC]
337 logical :: shelf_mass_is_dynamic, override_shelf_movement, active_shelf_dynamics
338 character(len=40) :: mdl = "MOM_ice_shelf_dyn" ! This module's name.
339 integer :: isd, ied, jsd, jed, isdb, iedb, jsdb, jedb
340
341 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
342 isdb = g%IsdB ; iedb = g%IedB ; jsdb = g%JsdB ; jedb = g%JedB
343
344 if (associated(cs)) then
345 call mom_error(fatal, "MOM_ice_shelf_dyn.F90, register_ice_shelf_dyn_restarts: "// &
346 "called with an associated control structure.")
347 return
348 endif
349 allocate(cs)
350
351 override_shelf_movement = .false. ; active_shelf_dynamics = .false.
352 call get_param(param_file, mdl, "DYNAMIC_SHELF_MASS", shelf_mass_is_dynamic, &
353 "If true, the ice sheet mass can evolve with time.", &
354 default=.false., do_not_log=.true.)
355 if (shelf_mass_is_dynamic) then
356 call get_param(param_file, mdl, "OVERRIDE_SHELF_MOVEMENT", override_shelf_movement, &
357 "If true, user provided code specifies the ice-shelf "//&
358 "movement instead of the dynamic ice model.", default=.false., do_not_log=.true.)
359 active_shelf_dynamics = .not.override_shelf_movement
360 endif
361
362 if (active_shelf_dynamics) then
363 call get_param(param_file, mdl, "MISSING_SHELF_TEMPERATURE", t_shelf_missing, &
364 "An ice shelf temperature to use where there is no ice shelf.",&
365 units="degC", default=-10.0, scale=us%degC_to_C, do_not_log=.true.)
366
367 call get_param(param_file, mdl, "NUMBER_OF_ICE_VISCOSITY_QUADRATURE_POINTS", cs%visc_qps, &
368 "Number of ice viscosity quadrature points. Either 1 (cell-centered) for 4", &
369 units="none", default=1)
370 if (cs%visc_qps/=1 .and. cs%visc_qps/=4) call mom_error (fatal, &
371 "NUMBER OF ICE_VISCOSITY_QUADRATURE_POINTS must be 1 or 4")
372
373 call get_param(param_file, mdl, "FIRST_DIRECTION_IS", cs%first_direction_IS, &
374 "An integer that indicates which direction goes first "//&
375 "in parts of the code that use directionally split "//&
376 "updates (e.g. advection), with even numbers (or 0) used for x- first "//&
377 "and odd numbers used for y-first.", default=0)
378 call get_param(param_file, mdl, "ALTERNATE_FIRST_DIRECTION_IS", cs%alternate_first_direction_IS, &
379 "If true, after every advection call, alternate whether the x- or y- "//&
380 "direction advection updates occur first. "//&
381 "If this is true, FIRST_DIRECTION applies at the start of a new run or if "//&
382 "the next first direction can not be found in the restart file.", default=.false.)
383
384 allocate(cs%u_shelf(isdb:iedb,jsdb:jedb), source=0.0)
385 allocate(cs%v_shelf(isdb:iedb,jsdb:jedb), source=0.0)
386 allocate(cs%t_shelf(isd:ied,jsd:jed), source=t_shelf_missing) ! [C ~> degC]
387 allocate(cs%ice_visc(isd:ied,jsd:jed,cs%visc_qps), source=0.0)
388 allocate(cs%newton_visc_factor(isd:ied,jsd:jed,cs%visc_qps), source=0.0)
389 allocate(cs%newton_str_ux(isd:ied,jsd:jed,cs%visc_qps), source=0.0)
390 allocate(cs%newton_str_vy(isd:ied,jsd:jed,cs%visc_qps), source=0.0)
391 allocate(cs%newton_str_sh(isd:ied,jsd:jed,cs%visc_qps), source=0.0)
392 allocate(cs%newton_umid(isd:ied,jsd:jed), source=0.0)
393 allocate(cs%newton_vmid(isd:ied,jsd:jed), source=0.0)
394 allocate(cs%newton_drag_coef(isd:ied,jsd:jed), source=0.0)
395 allocate(cs%AGlen_visc(isd:ied,jsd:jed), source=2.261e-25) ! [Pa-3 s-1]
396 allocate(cs%basal_traction(isd:ied,jsd:jed), source=0.0) ! [R Z L2 T-1 ~> kg s-1]
397 allocate(cs%C_basal_friction(isd:ied,jsd:jed), source=5.0e10*us%Pa_to_RLZ_T2)
398 ! Units of [R L Z T-2 (s m-1)^n_sliding ~> Pa (s m-1)^n_sliding]
399 allocate(cs%OD_av(isd:ied,jsd:jed), source=0.0)
400 allocate(cs%ground_frac(isd:ied,jsd:jed), source=0.0)
401 allocate(cs%taudx_shelf(isdb:iedb,jsdb:jedb), source=0.0)
402 allocate(cs%taudy_shelf(isdb:iedb,jsdb:jedb), source=0.0)
403 allocate(cs%sx_shelf(isd:ied,jsd:jed), source=0.0)
404 allocate(cs%sy_shelf(isd:ied,jsd:jed), source=0.0)
405 allocate(cs%bed_elev(isd:ied,jsd:jed), source=0.0)
406 allocate(cs%u_bdry_val(isdb:iedb,jsdb:jedb), source=0.0)
407 allocate(cs%v_bdry_val(isdb:iedb,jsdb:jedb), source=0.0)
408 allocate(cs%u_face_mask_bdry(isdb:iedb,jsdb:jedb), source=-2.0)
409 allocate(cs%v_face_mask_bdry(isdb:iedb,jsdb:jedb), source=-2.0)
410 allocate(cs%h_bdry_val(isd:ied,jsd:jed), source=0.0)
411
412 ! additional restarts for ice shelf state
413 call register_restart_field(cs%u_shelf, "u_shelf", .false., restart_cs, &
414 "ice sheet/shelf u-velocity", &
415 units="m s-1", conversion=us%L_T_to_m_s, hor_grid='Bu')
416 call register_restart_field(cs%v_shelf, "v_shelf", .false., restart_cs, &
417 "ice sheet/shelf v-velocity", &
418 units="m s-1", conversion=us%L_T_to_m_s, hor_grid='Bu')
419 call register_restart_field(cs%u_bdry_val, "u_bdry_val", .false., restart_cs, &
420 "ice sheet/shelf boundary u-velocity", &
421 units="m s-1", conversion=us%L_T_to_m_s, hor_grid='Bu')
422 call register_restart_field(cs%v_bdry_val, "v_bdry_val", .false., restart_cs, &
423 "ice sheet/shelf boundary v-velocity", &
424 units="m s-1", conversion=us%L_T_to_m_s, hor_grid='Bu')
425 call register_restart_field(cs%u_face_mask_bdry, "u_face_mask_bdry", .false., restart_cs, &
426 "ice sheet/shelf boundary u-mask", "nondim", hor_grid='Bu')
427 call register_restart_field(cs%v_face_mask_bdry, "v_face_mask_bdry", .false., restart_cs, &
428 "ice sheet/shelf boundary v-mask", "nondim", hor_grid='Bu')
429
430 call register_restart_field(cs%OD_av, "OD_av", .true., restart_cs, &
431 "Average open ocean depth in a cell", "m", conversion=us%Z_to_m)
432 call register_restart_field(cs%ground_frac, "ground_frac", .true., restart_cs, &
433 "fractional degree of grounding", "nondim")
434 call register_restart_field(cs%C_basal_friction, "C_basal_friction", .true., restart_cs, &
435 "basal sliding coefficients", "Pa (s m-1)^n_sliding", conversion=us%RLZ_T2_to_Pa)
436 call register_restart_field(cs%AGlen_visc, "AGlen_visc", .true., restart_cs, &
437 "ice-stiffness parameter", "Pa-3 s-1")
438 call register_restart_field(cs%h_bdry_val, "h_bdry_val", .false., restart_cs, &
439 "ice thickness at the boundary", "m", conversion=us%Z_to_m)
440 call register_restart_field(cs%bed_elev, "bed elevation", .true., restart_cs, &
441 "bed elevation", "m", conversion=us%Z_to_m)
442 call register_restart_field(cs%first_dir_restart_IS, "first_direction_IS", .false., restart_cs, &
443 "Indicator of the first direction in split ice shelf calculations.", "nondim")
444 endif
445
446end subroutine register_ice_shelf_dyn_restarts
447
448!> Initializes shelf model data, parameters and diagnostics
449subroutine initialize_ice_shelf_dyn(param_file, Time, ISS, CS, G, US, diag, new_sim, Cp_ice, &
450 Input_start_time, directory, solo_ice_sheet_in)
451 type(param_file_type), intent(in) :: param_file !< A structure to parse for run-time parameters
452 type(time_type), intent(inout) :: time !< The clock that that will indicate the model time
453 type(ice_shelf_state), intent(in) :: iss !< A structure with elements that describe
454 !! the ice-shelf state
455 type(ice_shelf_dyn_cs), pointer :: cs !< A pointer to the ice shelf dynamics control structure
456 type(ocean_grid_type), intent(inout) :: g !< The grid type describing the ice shelf grid.
457 type(unit_scale_type), intent(in) :: us !< A structure containing unit conversion factors
458 type(diag_ctrl), target, intent(in) :: diag !< A structure that is used to regulate the diagnostic output.
459 logical, intent(in) :: new_sim !< If true this is a new simulation, otherwise
460 !! has been started from a restart file.
461 real, intent(in) :: cp_ice !< Heat capacity of ice [Q C-1 ~> J kg-1 degC-1]
462 type(time_type), intent(in) :: input_start_time !< The start time of the simulation.
463 character(len=*), intent(in) :: directory !< The directory where the ice sheet energy file goes.
464 logical, optional, intent(in) :: solo_ice_sheet_in !< If present, this indicates whether
465 !! a solo ice-sheet driver.
466
467 ! Local variables
468 real :: t_shelf_bdry ! A default ice shelf temperature to use for ice flowing
469 ! in through open boundaries [C ~> degC]
470 !This include declares and sets the variable "version".
471# include "version_variable.h"
472 character(len=200) :: ic_file,filename,inputdir
473 character(len=40) :: var_name
474 character(len=40) :: mdl = "MOM_ice_shelf_dyn" ! This module's name.
475 logical :: shelf_mass_is_dynamic, override_shelf_movement, active_shelf_dynamics
476 logical :: enable_bugs ! If true, the defaults for recently added bug-fix flags are set to
477 ! recreate the bugs, or if false bugs are only used if actively selected.
478 logical :: debug
479 integer :: i, j, isd, ied, jsd, jed, isdq, iedq, jsdq, jedq, iters
480 character(len=200) :: is_energyfile ! The name of the energy file.
481 character(len=32) :: filename_appendix = '' ! FMS appendix to filename for ensemble runs
482
483 isdq = g%isdB ; iedq = g%iedB ; jsdq = g%jsdB ; jedq = g%jedB
484 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
485
486 if (.not.associated(cs)) then
487 call mom_error(fatal, "MOM_ice_shelf_dyn.F90, initialize_ice_shelf_dyn: "// &
488 "called with an associated control structure.")
489 return
490 endif
491 if (cs%module_is_initialized) then
492 call mom_error(warning, "MOM_ice_shelf_dyn.F90, initialize_ice_shelf_dyn was "//&
493 "called with a control structure that has already been initialized.")
494 endif
495 cs%module_is_initialized = .true.
496
497 cs%diag => diag ! ; CS%Time => Time
498
499 ! Read all relevant parameters and write them to the model log.
500 call log_version(param_file, mdl, version, "")
501 call get_param(param_file, mdl, "DEBUG", debug, default=.false.)
502 call get_param(param_file, mdl, "DEBUG_IS", cs%debug, &
503 "If true, write verbose debugging messages for the ice shelf.", &
504 default=debug)
505 call get_param(param_file, mdl, "DYNAMIC_SHELF_MASS", shelf_mass_is_dynamic, &
506 "If true, the ice sheet mass can evolve with time.", &
507 default=.false.)
508 override_shelf_movement = .false. ; active_shelf_dynamics = .false.
509 if (shelf_mass_is_dynamic) then
510 call get_param(param_file, mdl, "OVERRIDE_SHELF_MOVEMENT", override_shelf_movement, &
511 "If true, user provided code specifies the ice-shelf "//&
512 "movement instead of the dynamic ice model.", default=.false., do_not_log=.true.)
513 active_shelf_dynamics = .not.override_shelf_movement
514
515 call get_param(param_file, mdl, "GROUNDING_LINE_INTERPOLATE", cs%GL_regularize, &
516 "If true, regularize the floatation condition at the "//&
517 "grounding line as in Goldberg Holland Schoof 2009.", default=.false.)
518 call get_param(param_file, mdl, "GROUNDING_LINE_INTERP_SUBGRID_N", cs%n_sub_regularize, &
519 "The number of sub-partitions of each cell over which to "//&
520 "integrate for the interpolated grounding line. Each cell "//&
521 "is divided into NxN equally-sized rectangles, over which the "//&
522 "basal contribution is integrated by iterative quadrature.", &
523 default=0)
524 call get_param(param_file, mdl, "GROUNDING_LINE_COUPLE", cs%GL_couple, &
525 "If true, let the floatation condition be determined by "//&
526 "ocean column thickness. This means that update_OD_ffrac "//&
527 "will be called. GL_REGULARIZE and GL_COUPLE are exclusive.", &
528 default=.false., do_not_log=cs%GL_regularize)
529 if (cs%GL_regularize) cs%GL_couple = .false.
530 if (present(solo_ice_sheet_in)) then
531 if (solo_ice_sheet_in) cs%GL_couple = .false.
532 endif
533 if (cs%GL_regularize .and. (cs%n_sub_regularize == 0)) call mom_error (fatal, &
534 "GROUNDING_LINE_INTERP_SUBGRID_N must be a positive integer if GL regularization is used")
535 call get_param(param_file, mdl, "ICE_SHELF_CFL_FACTOR", cs%CFL_factor, &
536 "A factor used to limit timestep as CFL_FACTOR * min (\Delta x / u). "//&
537 "This is only used with an ice-only model.", units="nondim", default=0.25)
538 endif
539 call get_param(param_file, mdl, "RHO_0", cs%density_ocean_avg, &
540 "avg ocean density used in floatation cond", &
541 units="kg m-3", default=1035., scale=us%kg_m3_to_R)
542 if (active_shelf_dynamics) then
543 call get_param(param_file, mdl, "ICE_VELOCITY_TIMESTEP", cs%velocity_update_time_step, &
544 "seconds between ice velocity calcs", units="s", scale=us%s_to_T, &
545 fail_if_missing=.true.)
546 call get_param(param_file, mdl, "G_EARTH", cs%g_Earth, &
547 "The gravitational acceleration of the Earth.", &
548 units="m s-2", default=9.80, scale=us%m_s_to_L_T**2*us%Z_to_m)
549
550 call get_param(param_file, mdl, "MIN_H_SHELF", cs%min_h_shelf, &
551 "min. ice thickness used during ice dynamics", &
552 units="m", default=0.,scale=us%m_to_Z)
553 call get_param(param_file, mdl, "MIN_BASAL_TRACTION", cs%min_basal_traction, &
554 "min. allowed basal traction. Input is in [Pa m-1 yr], but is converted when read in to [Pa m-1 s]", &
555 units="Pa m-1 yr", default=0., scale=365.0*86400.0*us%Pa_to_RLZ_T2*us%L_T_to_m_s)
556 call get_param(param_file, mdl, "MAX_SURFACE_SLOPE", cs%max_surface_slope, &
557 "max. allowed ice-sheet surface slope. To ignore, set to zero.", &
558 units="none", default=0., scale=us%m_to_Z/us%m_to_L)
559 call get_param(param_file, mdl, "MIN_ICE_VISC", cs%min_ice_visc, &
560 "min. allowed Glen's law ice viscosity", &
561 units="Pa s", default=0., scale=us%Pa_to_RL2_T2*us%s_to_T)
562
563 call get_param(param_file, mdl, "GLEN_EXPONENT", cs%n_glen, &
564 "nonlinearity exponent in Glen's Law", &
565 units="none", default=3.)
566 call get_param(param_file, mdl, "MIN_STRAIN_RATE_GLEN", cs%eps_glen_min, &
567 "min. strain rate to avoid infinite Glen's law viscosity", &
568 units="s-1", default=1.e-19, scale=us%T_to_s)
569 call get_param(param_file, mdl, "BASAL_FRICTION_EXP", cs%n_basal_fric, &
570 "Exponent in sliding law \tau_b = C u^(n_basal_fric)", &
571 units="none", fail_if_missing=.true.)
572 call get_param(param_file, mdl, "USE_COULOMB_FRICTION", cs%CoulombFriction, &
573 "Use Coulomb Friction Law", &
574 units="none", default=.false., fail_if_missing=.false.)
575 call get_param(param_file, mdl, "CF_MinN", cs%CF_MinN, &
576 "Minimum Coulomb friction effective pressure", &
577 units="Pa", default=1.0, scale=us%Pa_to_RLZ_T2, fail_if_missing=.false.)
578 call get_param(param_file, mdl, "CF_PostPeak", cs%CF_PostPeak, &
579 "Coulomb friction post peak exponent", &
580 units="none", default=1.0, fail_if_missing=.false.)
581 call get_param(param_file, mdl, "CF_Max", cs%CF_Max, &
582 "Coulomb friction maximum coefficient", &
583 units="none", default=0.5, fail_if_missing=.false.)
584
585 call get_param(param_file, mdl, "DENSITY_ICE", cs%density_ice, &
586 "A typical density of ice.", units="kg m-3", default=917.0, scale=us%kg_m3_to_R)
587 call get_param(param_file, mdl, "CONJUGATE_GRADIENT_TOLERANCE", cs%cg_tolerance, &
588 "tolerance in CG solver, relative to initial residual", units="nondim", default=1.e-6)
589 cs%cg_tol_newton = cs%cg_tolerance ! Will be tightened adaptively during Newton iterations
590 call get_param(param_file, mdl, "ICE_NONLINEAR_TOLERANCE", cs%nonlinear_tolerance, &
591 "nonlin tolerance in iterative velocity solve", units="nondim", default=1.e-6)
592 call get_param(param_file, mdl, "NEWTON_AFTER_TOLERANCE", cs%newton_after_tolerance, &
593 "Switch from Picard to Newton iterations in the nonlinear ice velocity solve when "//&
594 "the fractional nonlinear residual falls below this tolerance.",&
595 units="none", default=cs%nonlinear_tolerance)
596 call get_param(param_file, mdl, "NEWTON_ADAPT_CG_TOL", cs%newton_adapt_cg_tol, &
597 "Use an adaptive CG tolerance during Newton iterations.", default=.true.)
598 call get_param(param_file, mdl, "CONJUGATE_GRADIENT_MAXIT", cs%cg_max_iterations, &
599 "max iteratiions in CG solver", default=2000)
600 call get_param(param_file, mdl, "THRESH_FLOAT_COL_DEPTH", cs%thresh_float_col_depth, &
601 "min ocean thickness to consider ice *floating*; "//&
602 "will only be important with use of tides", &
603 units="m", default=1.e-3, scale=us%m_to_Z)
604 call get_param(param_file, mdl, "NONLIN_SOLVE_ERR_MODE", cs%nonlin_solve_err_mode, &
605 "Choose whether nonlin error in vel solve is based on nonlinear "//&
606 "residual (1), relative change since last iteration (2), or change in norm (3)", default=3)
607
608 call get_param(param_file, mdl, "SHELF_MOVING_FRONT", cs%moving_shelf_front, &
609 "Specify whether to advance shelf front (and calve).", &
610 default=.false.)
611 call get_param(param_file, mdl, "CALVE_TO_MASK", cs%calve_to_mask, &
612 "If true, do not allow an ice shelf where prohibited by a mask.", &
613 default=.false.)
614 call get_param(param_file, mdl, "ADVECT_SHELF", cs%advect_shelf, &
615 "If true, advect ice shelf and evolve thickness", &
616 default=.true.)
617 call get_param(param_file, mdl, "REENTRANT_X", cs%reentrant_x, &
618 " If true, the domain is zonally reentrant.", &
619 default=.false.)
620 call get_param(param_file, mdl, "REENTRANT_Y", cs%reentrant_y, &
621 " If true, the domain is meridionally reentrant.", &
622 default=.false.)
623 call get_param(param_file, mdl, "ICE_VISCOSITY_COMPUTE", cs%ice_viscosity_compute, &
624 "If MODEL, compute ice viscosity internally using 1 or 4 quadrature points, "//&
625 "if OBS read from a file, "//&
626 "if CONSTANT a constant value (for debugging).", &
627 default="MODEL")
628
629 call get_param(param_file, mdl, "ENABLE_BUGS_BY_DEFAULT", enable_bugs, &
630 default=.true., do_not_log=.true.) ! This is logged from MOM.F90.
631 call get_param(param_file, mdl, "ICE_SHELF_TOP_SLOPE_BUG", cs%shelf_top_slope_bugs, &
632 "If true, use directionally inconsistent estimates of the grid spacing when "//&
633 "calculating the ice shelf surface slope, and underestimate slopes near the "//&
634 "edge of the ice shelf by a factor of 2.", default=enable_bugs)
635
636 if ((cs%visc_qps/=1) .and. (trim(cs%ice_viscosity_compute) /= "MODEL")) then
637 call mom_error(fatal, "NUMBER_OF_ICE_VISCOSITY_QUADRATURE_POINTS must be 1 unless ICE_VISCOSITY_COMPUTE==MODEL.")
638 endif
639 call get_param(param_file, mdl, "INFLOW_SHELF_TEMPERATURE", t_shelf_bdry, &
640 "A default ice shelf temperature to use for ice flowing in through "//&
641 "open boundaries.", units="degC", default=-15.0, scale=us%degC_to_C)
642 endif
643 call get_param(param_file, mdl, "MISSING_SHELF_TEMPERATURE", cs%T_shelf_missing, &
644 "An ice shelf temperature to use where there is no ice shelf.",&
645 units="degC", default=-10.0, scale=us%degC_to_C)
646 call get_param(param_file, mdl, "MIN_THICKNESS_SIMPLE_CALVE", cs%min_thickness_simple_calve, &
647 "Min thickness rule for the VERY simple calving law",&
648 units="m", default=0.0, scale=us%m_to_Z)
649 cs%Cp_ice = cp_ice !Heat capacity of ice (J kg-1 K-1), needed for heat flux of any bergs calved from
650 !the ice shelf and for ice sheet temperature solver
651 !for write_ice_shelf_energy
652 ! Note that the units of CS%Timeunit are the MKS units of [s].
653 call get_param(param_file, mdl, "TIMEUNIT", cs%Timeunit, &
654 "The time unit in seconds a number of input fields", &
655 units="s", default=86400.0)
656 if (cs%Timeunit < 0.0) cs%Timeunit = 86400.0
657 call get_param(param_file, mdl, "ENERGYSAVEDAYS",cs%energysavedays, &
658 "The interval in units of TIMEUNIT between saves of the "//&
659 "energies of the run and other globally summed diagnostics.",&
660 default=set_time(0,days=1), timeunit=cs%Timeunit)
661 call get_param(param_file, mdl, "ENERGYSAVEDAYS_GEOMETRIC",cs%energysavedays_geometric, &
662 "The starting interval in units of TIMEUNIT for the first call "//&
663 "to save the energies of the run and other globally summed diagnostics. "//&
664 "The interval increases by a factor of 2. after each call to write_ice_shelf_energy.",&
665 default=set_time(seconds=0), timeunit=cs%Timeunit)
666 if ((time_type_to_real(cs%energysavedays_geometric) > 0.) .and. &
667 (cs%energysavedays_geometric < cs%energysavedays)) then
668 cs%energysave_geometric = .true.
669 else
670 cs%energysave_geometric = .false.
671 endif
672 cs%Start_time = input_start_time
673 call get_param(param_file, mdl, "ICE_SHELF_ENERGYFILE", is_energyfile, &
674 "The file to use to write the energies and globally "//&
675 "summed diagnostics.", default="ice_shelf.stats")
676 !query fms_io if there is a filename_appendix (for ensemble runs)
677 call get_filename_appendix(filename_appendix)
678 if (len_trim(filename_appendix) > 0) then
679 is_energyfile = trim(is_energyfile) //'.'//trim(filename_appendix)
680 endif
681
682 cs%IS_energyfile = trim(slasher(directory))//trim(is_energyfile)
683 call log_param(param_file, mdl, "output_path/ENERGYFILE", cs%IS_energyfile)
684#ifdef STATSLABEL
685 cs%IS_energyfile = trim(cs%IS_energyfile)//"."//trim(adjustl(statslabel))
686#endif
687
688 ! Allocate memory in the ice shelf dynamics control structure that was not
689 ! previously allocated for registration for restarts.
690
691 if (active_shelf_dynamics) then
692 allocate( cs%t_bdry_val(isd:ied,jsd:jed), source=t_shelf_bdry) ! [C ~> degC]
693 allocate( cs%u_face_mask(isdq:iedq,jsdq:jedq), source=0.0)
694 allocate( cs%v_face_mask(isdq:iedq,jsdq:jedq), source=0.0)
695 allocate( cs%u_flux_bdry_val(isdq:iedq,jsd:jed), source=0.0)
696 allocate( cs%v_flux_bdry_val(isd:ied,jsdq:jedq), source=0.0)
697 allocate( cs%umask(isdq:iedq,jsdq:jedq), source=-1.0)
698 allocate( cs%vmask(isdq:iedq,jsdq:jedq), source=-1.0)
699 allocate( cs%tmask(isdq:iedq,jsdq:jedq), source=-1.0)
700 allocate( cs%float_cond(isd:ied,jsd:jed))
701
702 cs%OD_rt_counter = 0
703 allocate( cs%OD_rt(isd:ied,jsd:jed), source=0.0)
704 allocate( cs%ground_frac_rt(isd:ied,jsd:jed), source=0.0)
705
706 if (cs%calve_to_mask) then
707 allocate( cs%calve_mask(isd:ied,jsd:jed), source=0.0)
708 endif
709
710 allocate(cs%Phi(1:8,1:4,isd:ied,jsd:jed), source=0.0)
711 do j=g%jsd,g%jed ; do i=g%isd,g%ied
712 call bilinear_shape_fn_grid(g, i, j, cs%Phi(:,:,i,j))
713 enddo ; enddo
714
715 if (cs%GL_regularize) then
716 allocate(cs%Phisub(2,2,cs%n_sub_regularize,cs%n_sub_regularize,2,2), source=0.0)
717 call bilinear_shape_functions_subgrid(cs%Phisub, cs%n_sub_regularize)
718 endif
719
720 if ((trim(cs%ice_viscosity_compute) == "MODEL") .and. cs%visc_qps==1) then
721 !for calculating viscosity and 1 cell-centered quadrature point per cell
722 allocate(cs%PhiC(1:8,g%isc:g%iec,g%jsc:g%jec), source=0.0)
723 do j=g%jsc,g%jec ; do i=g%isc,g%iec
724 call bilinear_shape_fn_grid_1qp(g, i, j, cs%PhiC(:,i,j))
725 enddo ; enddo
726 endif
727
728 cs%elapsed_velocity_time = 0.0
729
730 call update_velocity_masks(cs, g, iss%hmask, cs%umask, cs%vmask, cs%u_face_mask, cs%v_face_mask)
731 endif
732
733 ! Take additional initialization steps, for example of dependent variables.
734 if (active_shelf_dynamics .and. .not.new_sim) then
735
736 call pass_var(cs%OD_av,g%domain, complete=.false.)
737 call pass_var(cs%ground_frac, g%domain, complete=.false.)
738 call pass_var(cs%basal_traction, g%domain, complete=.false.)
739 call pass_var(cs%AGlen_visc, g%domain, complete=.false.)
740 call pass_var(cs%bed_elev, g%domain, complete=.false.)
741 call pass_var(cs%C_basal_friction, g%domain, complete=.false.)
742 call pass_var(cs%h_bdry_val, g%domain, complete=.true.)
743 call pass_var(cs%ice_visc, g%domain)
744
745 call pass_vector(cs%u_bdry_val, cs%v_bdry_val, g%domain, to_all, bgrid_ne, complete=.false.)
746 call pass_vector(cs%u_face_mask_bdry, cs%v_face_mask_bdry, g%domain, to_all, bgrid_ne, complete=.true.)
747 call update_velocity_masks(cs, g, iss%hmask, cs%umask, cs%vmask, cs%u_face_mask, cs%v_face_mask)
748
749 ! This is unfortunately necessary (?); if grid is not symmetric the boundary values
750 ! of u and v are otherwise not set till the end of the first linear solve, and so
751 ! viscosity is not calculated correctly.
752 ! This has to occur after init_boundary_values or some of the arrays on the
753 ! right hand side have not been set up yet.
754 if (.not. g%symmetric) then
755 do j=g%jsd,g%jed ; do i=g%isd,g%ied
756 if ((i+g%idg_offset) == (g%domain%nihalo+1)) then
757 if (cs%u_face_mask(i-1,j) == 3) then
758 cs%u_shelf(i-1,j-1) = cs%u_bdry_val(i-1,j-1)
759 cs%u_shelf(i-1,j) = cs%u_bdry_val(i-1,j)
760 cs%v_shelf(i-1,j-1) = cs%v_bdry_val(i-1,j-1)
761 cs%v_shelf(i-1,j) = cs%v_bdry_val(i-1,j)
762 elseif (cs%u_face_mask(i-1,j) == 5) then
763 cs%u_shelf(i-1,j-1) = cs%u_bdry_val(i-1,j-1)
764 cs%u_shelf(i-1,j) = cs%u_bdry_val(i-1,j)
765 elseif (cs%u_face_mask(i-1,j) == 6) then
766 cs%v_shelf(i-1,j-1) = cs%v_bdry_val(i-1,j-1)
767 cs%v_shelf(i-1,j) = cs%v_bdry_val(i-1,j)
768 endif
769 endif
770 if ((j+g%jdg_offset) == (g%domain%njhalo+1)) then
771 if (cs%v_face_mask(i,j-1) == 3) then
772 cs%v_shelf(i-1,j-1) = cs%v_bdry_val(i-1,j-1)
773 cs%v_shelf(i,j-1) = cs%v_bdry_val(i,j-1)
774 cs%u_shelf(i-1,j-1) = cs%u_bdry_val(i-1,j-1)
775 cs%u_shelf(i,j-1) = cs%u_bdry_val(i,j-1)
776 elseif (cs%v_face_mask(i,j-1) == 5) then
777 cs%v_shelf(i-1,j-1) = cs%v_bdry_val(i-1,j-1)
778 cs%v_shelf(i,j-1) = cs%v_bdry_val(i,j-1)
779 elseif (cs%v_face_mask(i,j-1) == 6) then
780 cs%u_shelf(i-1,j-1) = cs%u_bdry_val(i-1,j-1)
781 cs%u_shelf(i,j-1) = cs%u_bdry_val(i,j-1)
782 endif
783 endif
784 enddo ; enddo
785 endif
786 call pass_vector(cs%u_shelf, cs%v_shelf, g%domain, to_all, bgrid_ne)
787 endif
788
789 if (active_shelf_dynamics) then
790 if (cs%first_dir_restart_IS > -1.0) then
791 cs%first_direction_IS = modulo(nint(cs%first_dir_restart_IS), 2)
792 else
793 cs%first_dir_restart_IS = real(modulo(cs%first_direction_IS, 2))
794 endif
795
796 ! If we are calving to a mask, i.e. if a mask exists where a shelf cannot, read the mask from a file.
797 if (cs%calve_to_mask) then
798 call mom_mesg(" MOM_ice_shelf.F90, initialize_ice_shelf: reading calving_mask")
799
800 call get_param(param_file, mdl, "INPUTDIR", inputdir, default=".")
801 inputdir = slasher(inputdir)
802 call get_param(param_file, mdl, "CALVING_MASK_FILE", ic_file, &
803 "The file with a mask for where calving might occur.", &
804 default="ice_shelf_h.nc")
805 call get_param(param_file, mdl, "CALVING_MASK_VARNAME", var_name, &
806 "The variable to use in masking calving.", &
807 default="area_shelf_h")
808
809 filename = trim(inputdir)//trim(ic_file)
810 call log_param(param_file, mdl, "INPUTDIR/CALVING_MASK_FILE", filename)
811 if (.not.file_exists(filename, g%Domain)) call mom_error(fatal, &
812 " calving mask file: Unable to open "//trim(filename))
813
814 call mom_read_data(filename,trim(var_name),cs%calve_mask,g%Domain)
815 do j=g%jsc,g%jec ; do i=g%isc,g%iec
816 if (cs%calve_mask(i,j) > 0.0) cs%calve_mask(i,j) = 1.0
817 enddo ; enddo
818 call pass_var(cs%calve_mask,g%domain)
819 endif
820
821 ! initialize basal friction coefficients
822 if (new_sim) then
823 call initialize_ice_c_basal_friction(cs%C_basal_friction, g, us, param_file)
824 call pass_var(cs%C_basal_friction, g%domain, complete=.false.)
825
826 ! initialize ice-stiffness AGlen
827 call initialize_ice_aglen(cs%AGlen_visc, cs%ice_viscosity_compute, g, us, param_file)
828 call pass_var(cs%AGlen_visc, g%domain, complete=.false.)
829
830 !initialize boundary conditions
831 call initialize_ice_shelf_boundary_from_file(cs%u_face_mask_bdry, cs%v_face_mask_bdry, &
832 cs%u_bdry_val, cs%v_bdry_val, cs%umask, cs%vmask, cs%h_bdry_val, &
833 iss%hmask, iss%h_shelf, g, us, param_file )
834 call pass_var(iss%hmask, g%domain, complete=.false.)
835 call pass_var(cs%h_bdry_val, g%domain, complete=.true.)
836 call pass_vector(cs%u_bdry_val, cs%v_bdry_val, g%domain, to_all, bgrid_ne, complete=.false.)
837 call pass_vector(cs%u_face_mask_bdry, cs%v_face_mask_bdry, g%domain, to_all, bgrid_ne, complete=.false.)
838
839 !initialize ice flow characteristic (velocities, bed elevation under the grounded part, etc) from file
840 call initialize_ice_flow_from_file(cs%bed_elev,cs%u_shelf, cs%v_shelf, cs%ground_frac, &
841 g, us, param_file)
842 call pass_vector(cs%u_shelf, cs%v_shelf, g%domain, to_all, bgrid_ne, complete=.true.)
843 call pass_var(cs%ground_frac, g%domain, complete=.false.)
844 call pass_var(cs%bed_elev, g%domain, complete=.true.)
845 call update_velocity_masks(cs, g, iss%hmask, cs%umask, cs%vmask, cs%u_face_mask, cs%v_face_mask)
846
847 do j=jsdq,jedq ; do i=isdq,iedq
848 if (cs%umask(i,j) == 3) then
849 cs%u_shelf(i,j) = cs%u_bdry_val(i,j)
850 elseif (cs%umask(i,j) == 0) then
851 cs%u_shelf(i,j) = 0
852 endif
853 if (cs%vmask(i,j) == 3) then
854 cs%v_shelf(i,j) = cs%v_bdry_val(i,j)
855 elseif (cs%vmask(i,j) == 0) then
856 cs%v_shelf(i,j) = 0
857 endif
858 enddo ; enddo
859 endif
860
861 ! Register diagnostics.
862 cs%id_u_shelf = register_diag_field('ice_shelf_model','u_shelf',cs%diag%axesB1, time, &
863 'x-velocity of ice', 'm yr-1', conversion=365.0*86400.0*us%L_T_to_m_s)
864 cs%id_v_shelf = register_diag_field('ice_shelf_model','v_shelf',cs%diag%axesB1, time, &
865 'y-velocity of ice', 'm yr-1', conversion=365.0*86400.0*us%L_T_to_m_s)
866 cs%id_shelf_speed = register_diag_field('ice_shelf_model','shelf_speed',cs%diag%axesB1, time, &
867 'speed of of ice shelf', 'm yr-1', conversion=365.0*86400.0*us%L_T_to_m_s)
868 cs%id_taudx_shelf = register_diag_field('ice_shelf_model','taudx_shelf',cs%diag%axesB1, time, &
869 'x-driving stress of ice', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
870 cs%id_taudy_shelf = register_diag_field('ice_shelf_model','taudy_shelf',cs%diag%axesB1, time, &
871 'y-driving stress of ice', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
872 cs%id_taud_shelf = register_diag_field('ice_shelf_model','taud_shelf',cs%diag%axesB1, time, &
873 'magnitude of driving stress of ice', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
874 cs%id_sx_shelf = register_diag_field('ice_shelf_model', 'sx_shelf', cs%diag%axesT1, time, &
875 'x-surface slope of ice', 'none')
876 cs%id_sy_shelf = register_diag_field('ice_shelf_model', 'sy_shelf', cs%diag%axesT1, time, &
877 'y-surface slope of ice', 'none')
878 cs%id_surf_slope_mag_shelf = register_diag_field('ice_shelf_model', 'surf_slope_mag_shelf', cs%diag%axesT1, time, &
879 'magnitude of surface slope of ice', 'none')
880 cs%id_u_mask = register_diag_field('ice_shelf_model','u_mask',cs%diag%axesB1, time, &
881 'mask for u-nodes', 'none')
882 cs%id_v_mask = register_diag_field('ice_shelf_model','v_mask',cs%diag%axesB1, time, &
883 'mask for v-nodes', 'none')
884 cs%id_ground_frac = register_diag_field('ice_shelf_model','ice_ground_frac',cs%diag%axesT1, time, &
885 'fraction of cell that is grounded', 'none')
886 cs%id_float_cond = register_diag_field('ice_shelf_model','float_cond',cs%diag%axesT1, time, &
887 'sub-cell grounding cells', 'none')
888 cs%id_col_thick = register_diag_field('ice_shelf_model','col_thick',cs%diag%axesT1, time, &
889 'ocean column thickness passed to ice model', 'm', conversion=us%Z_to_m)
890 cs%id_visc_shelf = register_diag_field('ice_shelf_model','ice_visc',cs%diag%axesT1, time, &
891 'vi-viscosity', 'Pa m s', conversion=us%RL2_T2_to_Pa*us%Z_to_m*us%T_to_s) !vertically integrated viscosity
892 cs%id_taub = register_diag_field('ice_shelf_model','taub_beta',cs%diag%axesT1, time, &
893 'taub', units='MPa yr m-1', conversion=1e-6*us%RLZ_T2_to_Pa/(365.0*86400.0*us%L_T_to_m_s))
894 cs%id_OD_av = register_diag_field('ice_shelf_model','OD_av',cs%diag%axesT1, time, &
895 'intermediate ocean column thickness passed to ice model', 'm', conversion=us%Z_to_m)
896
897 cs%id_duHdx = register_diag_field('ice_shelf_model','duHdx',cs%diag%axesT1, time, &
898 'x-component of ice-sheet flux divergence', 'm yr-1', conversion=365.0*86400.0*us%Z_to_m*us%s_to_T)
899 cs%id_dvHdy = register_diag_field('ice_shelf_model','dvHdy',cs%diag%axesT1, time, &
900 'y-component of ice-sheet flux divergence', 'm yr-1', conversion=365.0*86400.0*us%Z_to_m*us%s_to_T)
901 cs%id_fluxdiv = register_diag_field('ice_shelf_model','fluxdiv',cs%diag%axesT1, time, &
902 'ice-sheet flux divergence', 'm yr-1', conversion=365.0*86400.0*us%Z_to_m*us%s_to_T)
903 cs%id_strainrate_xx = register_diag_field('ice_shelf_model','strainrate_xx',cs%diag%axesT1, time, &
904 'x-component of ice-shelf strain-rate', 'yr-1', conversion=365.0*86400.0*us%s_to_T)
905 cs%id_strainrate_yy = register_diag_field('ice_shelf_model','strainrate_yy',cs%diag%axesT1, time, &
906 'y-component of ice-shelf strain-rate', 'yr-1', conversion=365.0*86400.0*us%s_to_T)
907 cs%id_strainrate_xy = register_diag_field('ice_shelf_model','strainrate_xy',cs%diag%axesT1, time, &
908 'xy-component of ice-shelf strain-rate', 'yr-1', conversion=365.0*86400.0*us%s_to_T)
909 cs%id_pstrainrate_1 = register_diag_field('ice_shelf_model','pstrainrate_1',cs%diag%axesT1, time, &
910 'max principal horizontal ice-shelf strain-rate', 'yr-1', conversion=365.0*86400.0*us%s_to_T)
911 cs%id_pstrainrate_2 = register_diag_field('ice_shelf_model','pstrainrate_2',cs%diag%axesT1, time, &
912 'min principal horizontal ice-shelf strain-rate', 'yr-1', conversion=365.0*86400.0*us%s_to_T)
913 cs%id_devstress_xx = register_diag_field('ice_shelf_model','devstress_xx',cs%diag%axesT1, time, &
914 'x-component of ice-shelf deviatoric stress', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
915 cs%id_devstress_yy = register_diag_field('ice_shelf_model','devstress_yy',cs%diag%axesT1, time, &
916 'y-component of ice-shelf deviatoric stress', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
917 cs%id_devstress_xy = register_diag_field('ice_shelf_model','devstress_xy',cs%diag%axesT1, time, &
918 'xy-component of ice-shelf deviatoric stress', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
919 cs%id_pdevstress_1 = register_diag_field('ice_shelf_model','pdevstress_1',cs%diag%axesT1, time, &
920 'max principal horizontal ice-shelf deviatoric stress', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
921 cs%id_pdevstress_2 = register_diag_field('ice_shelf_model','pdevstress_2',cs%diag%axesT1, time, &
922 'min principal ice-shelf deviatoric stress', 'kPa', conversion=1.e-3*us%RLZ_T2_to_Pa)
923
924 !Update these variables so that they are nonzero in case
925 !IS_dynamics_post_data is called before update_ice_shelf
926 if (cs%id_taudx_shelf>0 .or. cs%id_taudy_shelf>0) &
927 call calc_shelf_driving_stress(cs, iss, g, us, cs%taudx_shelf, cs%taudy_shelf, cs%OD_av)
928 if (cs%id_taub>0) &
929 call calc_shelf_taub(cs, iss, g, us, cs%u_shelf, cs%v_shelf)
930 if (cs%id_visc_shelf>0) &
931 call calc_shelf_visc(cs, iss, g, us, cs%u_shelf, cs%v_shelf)
932 endif
933
934 if (new_sim) then
935 call update_od_ffrac_uncoupled(cs, g, iss%h_shelf(:,:))
936 endif
937
938end subroutine initialize_ice_shelf_dyn
939
940
941subroutine initialize_diagnostic_fields(CS, ISS, G, US, Time)
942 type(ice_shelf_dyn_cs), intent(inout) :: CS !< A pointer to the ice shelf control structure
943 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
944 !! the ice-shelf state
945 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
946 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
947 type(time_type), intent(in) :: Time !< The current model time
948
949 integer :: i, j, iters, isd, ied, jsd, jed
950 real :: rhoi_rhow
951 real :: OD ! Depth of open water below the ice shelf [Z ~> m]
952 type(time_type) :: dummy_time
953!
954 rhoi_rhow = cs%density_ice / cs%density_ocean_avg
955 dummy_time = set_time(0,0)
956 isd=g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
957
958 do j=jsd,jed
959 do i=isd,ied
960 od = cs%bed_elev(i,j) - rhoi_rhow * max(iss%h_shelf(i,j),cs%min_h_shelf)
961 if (od >= 0) then
962 ! ice thickness does not take up whole ocean column -> floating
963 cs%OD_av(i,j) = od
964 cs%ground_frac(i,j) = 0.
965 else
966 cs%OD_av(i,j) = 0.
967 cs%ground_frac(i,j) = 1.
968 endif
969 enddo
970 enddo
971
972 call ice_shelf_solve_outer(cs, iss, g, us, cs%u_shelf, cs%v_shelf,cs%taudx_shelf,cs%taudy_shelf, iters, time)
973end subroutine initialize_diagnostic_fields
974
975!> This function returns the global maximum advective timestep that can be taken based on the current
976!! ice velocities. Because it involves finding a global minimum, it can be surprisingly expensive.
977function ice_time_step_cfl(CS, ISS, G)
978 type(ice_shelf_dyn_cs), intent(inout) :: cs !< The ice shelf dynamics control structure
979 type(ice_shelf_state), intent(inout) :: iss !< A structure with elements that describe
980 !! the ice-shelf state
981 type(ocean_grid_type), intent(inout) :: g !< The grid structure used by the ice shelf.
982 real :: ice_time_step_cfl !< The maximum permitted timestep based on the ice velocities [T ~> s].
983
984 real :: dt_local, min_dt ! These should be the minimum stable timesteps at a CFL of 1 [T ~> s]
985 real :: min_vel ! A minimal velocity for estimating a timestep [L T-1 ~> m s-1]
986 integer :: i, j
987
988 min_dt = 5.0e17*g%US%s_to_T ! The starting maximum is roughly the lifetime of the universe.
989 min_vel = (1.0e-12/(365.0*86400.0)) * g%US%m_s_to_L_T
990 do j=g%jsc,g%jec ; do i=g%isc,g%iec ; if (iss%hmask(i,j) == 1.0 .or. iss%hmask(i,j)==3) then
991 dt_local = 2.0*g%areaT(i,j) / &
992 (((g%dyCu(i,j) * max(abs(cs%u_shelf(i,j) + cs%u_shelf(i,j-1)), min_vel)) + &
993 (g%dyCu(i-1,j)* max(abs(cs%u_shelf(i-1,j)+ cs%u_shelf(i-1,j-1)), min_vel))) + &
994 ((g%dxCv(i,j) * max(abs(cs%v_shelf(i,j) + cs%v_shelf(i-1,j)), min_vel)) + &
995 (g%dxCv(i,j-1)* max(abs(cs%v_shelf(i,j-1)+ cs%v_shelf(i-1,j-1)), min_vel))))
996
997 min_dt = min(min_dt, dt_local)
998 endif ; enddo ; enddo ! i- and j- loops
999
1000 call min_across_pes(min_dt)
1001
1002 ice_time_step_cfl = cs%CFL_factor * min_dt
1003
1004end function ice_time_step_cfl
1005
1006!> This subroutine updates the ice shelf velocities, mass, stresses and properties due to the
1007!! ice shelf dynamics.
1008subroutine update_ice_shelf(CS, ISS, G, US, time_step, Time, calve_ice_shelf_bergs, &
1009 ocean_mass, coupled_grounding, must_update_vel)
1010 type(ice_shelf_dyn_cs), intent(inout) :: cs !< The ice shelf dynamics control structure
1011 type(ice_shelf_state), intent(inout) :: iss !< A structure with elements that describe
1012 !! the ice-shelf state
1013 type(ocean_grid_type), intent(inout) :: g !< The grid structure used by the ice shelf.
1014 type(unit_scale_type), intent(in) :: us !< A structure containing unit conversion factors
1015 real, intent(in) :: time_step !< time step [T ~> s]
1016 type(time_type), intent(in) :: time !< The current model time
1017 logical, intent(in) :: calve_ice_shelf_bergs !< To convert ice flux through front
1018 !! to bergs
1019 real, dimension(SZDI_(G),SZDJ_(G)), &
1020 optional, intent(in) :: ocean_mass !< If present this is the mass per unit area
1021 !! of the ocean [R Z ~> kg m-2].
1022 logical, optional, intent(in) :: coupled_grounding !< If true, the grounding line is
1023 !! determined by coupled ice-ocean dynamics
1024 logical, optional, intent(in) :: must_update_vel !< Always update the ice velocities if true.
1025 integer :: iters
1026 logical :: update_ice_vel, coupled_gl
1027
1028 update_ice_vel = .false.
1029 if (present(must_update_vel)) update_ice_vel = must_update_vel
1030
1031 coupled_gl = .false.
1032 if (present(ocean_mass) .and. present(coupled_grounding)) coupled_gl = coupled_grounding
1033!
1034 if (cs%advect_shelf) then
1035 call ice_shelf_advect(cs, iss, g, time_step, time, calve_ice_shelf_bergs)
1036 if (cs%alternate_first_direction_IS) then
1037 cs%first_direction_IS = modulo(cs%first_direction_IS+1,2)
1038 cs%first_dir_restart_IS = real(cs%first_direction_IS)
1039 endif
1040 endif
1041 cs%elapsed_velocity_time = cs%elapsed_velocity_time + time_step
1042 if (cs%elapsed_velocity_time >= cs%velocity_update_time_step) update_ice_vel = .true.
1043
1044 if (coupled_gl) then
1045 call update_od_ffrac(cs, g, us, ocean_mass, update_ice_vel)
1046 elseif (update_ice_vel) then
1047 call update_od_ffrac_uncoupled(cs, g, iss%h_shelf(:,:))
1048 cs%GL_couple=.false.
1049 endif
1050
1051 if (update_ice_vel) then
1052 call ice_shelf_solve_outer(cs, iss, g, us, cs%u_shelf, cs%v_shelf,cs%taudx_shelf,cs%taudy_shelf, iters, time)
1053 cs%elapsed_velocity_time = 0.0
1054 endif
1055
1056! call ice_shelf_temp(CS, ISS, G, US, time_step, ISS%water_flux, Time)
1057
1058end subroutine update_ice_shelf
1059
1060subroutine volume_above_floatation(CS, G, ISS, vaf, hemisphere)
1061 type(ice_shelf_dyn_cs), intent(in) :: cs !< The ice shelf dynamics control structure
1062 type(ocean_grid_type), intent(in) :: g !< The grid structure used by the ice shelf.
1063 type(ice_shelf_state), intent(in) :: iss !< A structure with elements that describe
1064 !! the ice-shelf state
1065 real, intent(out) :: vaf !< area integrated volume above floatation [Z L2 ~> m3]
1066 integer, optional, intent(in) :: hemisphere !< 0 for Antarctica only, 1 for Greenland only. Otherwise, all ice sheets
1067 integer :: is_id ! local copy of hemisphere
1068 real, dimension(SZI_(G),SZJ_(G)) :: vaf_cell !< cell-wise volume above floatation [Z L2 ~> m3]
1069 integer, dimension(SZI_(G),SZJ_(G)) :: mask ! a mask for active cells depending on hemisphere indicated
1070 integer :: is,ie,js,je,i,j
1071 real :: rhoi_rhow, rhow_rhoi
1072
1073 if (cs%GL_couple) &
1074 call mom_error(fatal, "MOM_ice_shelf_dyn, volume above floatation calculation assumes GL_couple=.FALSE..")
1075
1076 rhoi_rhow = cs%density_ice / cs%density_ocean_avg
1077 rhow_rhoi = cs%density_ocean_avg / cs%density_ice
1078 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1079
1080 if (present(hemisphere)) then
1081 is_id=hemisphere
1082 else
1083 is_id=-1
1084 endif
1085
1086 mask(:,:)=0
1087 if (is_id==0) then !Antarctica (S. Hemisphere) only
1088 do j = js,je ; do i = is,ie
1089 if (iss%hmask(i,j)>0 .and. g%geoLatT(i,j)<=0.0) mask(i,j)=1
1090 enddo ; enddo
1091 elseif (is_id==1) then !Greenland (N. Hemisphere) only
1092 do j = js,je ; do i = is,ie
1093 if (iss%hmask(i,j)>0 .and. g%geoLatT(i,j)>0.0) mask(i,j)=1
1094 enddo ; enddo
1095 else !All ice sheets
1096 mask(is:ie,js:je)=iss%hmask(is:ie,js:je)
1097 endif
1098
1099 vaf_cell(:,:)=0.0
1100 do j = js,je ; do i = is,ie
1101 if (mask(i,j)>0) then
1102 if (cs%bed_elev(i,j) <= 0) then
1103 !grounded above sea level
1104 vaf_cell(i,j) = iss%h_shelf(i,j) * iss%area_shelf_h(i,j)
1105 else
1106 !grounded if vaf_cell(i,j) > 0
1107 vaf_cell(i,j) = max(iss%h_shelf(i,j) - rhow_rhoi * cs%bed_elev(i,j), 0.0) * iss%area_shelf_h(i,j)
1108 endif
1109 endif
1110 enddo ; enddo
1111
1112 vaf = reproducing_sum(vaf_cell, unscale=g%US%Z_to_m*g%US%L_to_m**2)
1113end subroutine volume_above_floatation
1114
1115!> multiplies a variable with the ice sheet grounding fraction
1116subroutine masked_var_grounded(G,CS,var,varout)
1117 type(ocean_grid_type), intent(in) :: g !< The grid structure used by the ice shelf.
1118 type(ice_shelf_dyn_cs), intent(in) :: cs !< The ice shelf dynamics control structure
1119 real, dimension(SZI_(G),SZJ_(G)), intent(in) :: var !< variable in
1120 real, dimension(SZI_(G),SZJ_(G)), intent(out) :: varout !<variable out
1121 integer :: i, j
1122 do j = g%jsc,g%jec ; do i = g%isc,g%iec
1123 varout(i,j) = var(i,j) * cs%ground_frac(i,j)
1124 enddo ; enddo
1125end subroutine masked_var_grounded
1126
1127!> Ice shelf dynamics post_data calls
1128subroutine is_dynamics_post_data(time_step, Time, CS, ISS, G)
1129 real :: time_step !< Length of time for post data averaging [T ~> s].
1130 type(time_type), intent(in) :: time !< The current model time
1131 type(ice_shelf_dyn_cs), intent(inout) :: cs !< The ice shelf dynamics control structure
1132 type(ice_shelf_state), intent(inout) :: iss !< A structure with elements that describe
1133 !! the ice-shelf state
1134 type(ocean_grid_type), intent(in) :: g !< The grid structure used by the ice shelf.
1135 real, dimension(SZDIB_(G),SZDJB_(G)) :: taud_x, taud_y, taud ! area-averaged driving stress [R L2 T-2 ~> Pa]
1136 real, dimension(SZDI_(G),SZDJ_(G)) :: ice_visc ! area-averaged vertically integrated ice viscosity
1137 !! [R L2 Z T-1 ~> Pa s m]
1138 real, dimension(SZDI_(G),SZDJ_(G)) :: basal_tr ! area-averaged taub_beta field related to basal traction,
1139 !! [R L T-1 ~> Pa s m-1]
1140 real, dimension(SZDI_(G),SZDJ_(G)) :: surf_slope ! the surface slope of the ice shelf/sheet [nondim]
1141 real, dimension(SZDIB_(G),SZDJB_(G)) :: ice_speed ! ice sheet flow speed [L T-1 ~> m s-1]
1142
1143 integer :: i, j
1144
1145 call enable_averages(time_step, time, cs%diag)
1146 if (cs%id_col_thick > 0) call post_data(cs%id_col_thick, cs%OD_av, cs%diag)
1147 if (cs%id_u_shelf > 0) call post_data(cs%id_u_shelf, cs%u_shelf, cs%diag)
1148 if (cs%id_v_shelf > 0) call post_data(cs%id_v_shelf, cs%v_shelf, cs%diag)
1149 if (cs%id_shelf_speed > 0) then
1150 do j=g%jscB,g%jecB ; do i=g%iscB,g%iecB
1151 ice_speed(i,j) = sqrt((cs%u_shelf(i,j)**2) + (cs%v_shelf(i,j)**2))
1152 enddo ; enddo
1153 call post_data(cs%id_shelf_speed, ice_speed, cs%diag)
1154 endif
1155! if (CS%id_t_shelf > 0) call post_data(CS%id_t_shelf, CS%t_shelf, CS%diag)
1156 if (cs%id_taudx_shelf > 0) then
1157 do j=g%jscB,g%jecB ; do i=g%iscB,g%iecB
1158 taud_x(i,j) = cs%taudx_shelf(i,j)*g%IareaBu(i,j)
1159 enddo ; enddo
1160 call post_data(cs%id_taudx_shelf, taud_x, cs%diag)
1161 endif
1162 if (cs%id_taudy_shelf > 0) then
1163 do j=g%jscB,g%jecB ; do i=g%iscB,g%iecB
1164 taud_y(i,j) = cs%taudy_shelf(i,j)*g%IareaBu(i,j)
1165 enddo ; enddo
1166 call post_data(cs%id_taudy_shelf, taud_y, cs%diag)
1167 endif
1168 if (cs%id_taud_shelf > 0) then
1169 do j=g%jscB,g%jecB ; do i=g%iscB,g%iecB
1170 taud(i,j) = sqrt((cs%taudx_shelf(i,j)**2)+(cs%taudy_shelf(i,j)**2))*g%IareaBu(i,j)
1171 enddo ; enddo
1172 call post_data(cs%id_taud_shelf, taud, cs%diag)
1173 endif
1174 if (cs%id_sx_shelf > 0) call post_data(cs%id_sx_shelf, cs%sx_shelf, cs%diag)
1175 if (cs%id_sy_shelf > 0) call post_data(cs%id_sy_shelf, cs%sy_shelf, cs%diag)
1176 if (cs%id_surf_slope_mag_shelf > 0) then
1177 do j=g%jsc,g%jec ; do i=g%isc,g%iec
1178 surf_slope(i,j) = sqrt((cs%sx_shelf(i,j)**2)+(cs%sy_shelf(i,j)**2))
1179 enddo ; enddo
1180 call post_data(cs%id_surf_slope_mag_shelf, surf_slope, cs%diag)
1181 endif
1182 if (cs%id_ground_frac > 0) call post_data(cs%id_ground_frac, cs%ground_frac, cs%diag)
1183 if (cs%id_float_cond > 0) call post_data(cs%id_float_cond, cs%float_cond, cs%diag)
1184 if (cs%id_OD_av >0) call post_data(cs%id_OD_av, cs%OD_av,cs%diag)
1185 if (cs%id_visc_shelf > 0) then
1186 call ice_visc_diag(cs,g,ice_visc)
1187 call post_data(cs%id_visc_shelf, ice_visc, cs%diag)
1188 endif
1189 if (cs%id_taub > 0) then
1190 do j=g%jsc,g%jec ; do i=g%isc,g%iec
1191 basal_tr(i,j) = cs%basal_traction(i,j)*g%IareaT(i,j)
1192 enddo ; enddo
1193 call post_data(cs%id_taub, basal_tr, cs%diag)
1194 endif
1195 if (cs%id_u_mask > 0) call post_data(cs%id_u_mask, cs%umask, cs%diag)
1196 if (cs%id_v_mask > 0) call post_data(cs%id_v_mask, cs%vmask, cs%diag)
1197 if (cs%id_ufb_mask > 0) call post_data(cs%id_ufb_mask, cs%u_face_mask_bdry, cs%diag)
1198 if (cs%id_vfb_mask > 0) call post_data(cs%id_vfb_mask, cs%v_face_mask_bdry, cs%diag)
1199! if (CS%id_t_mask > 0) call post_data(CS%id_t_mask, CS%tmask, CS%diag)
1200
1201 if (cs%id_duHdx > 0 .or. cs%id_dvHdy > 0 .or. cs%id_fluxdiv > 0 .or. &
1202 cs%id_devstress_xx > 0 .or. cs%id_devstress_yy > 0 .or. cs%id_devstress_xy > 0 .or. &
1203 cs%id_strainrate_xx > 0 .or. cs%id_strainrate_yy > 0 .or. cs%id_strainrate_xy > 0 .or. &
1204 cs%id_pdevstress_1 > 0 .or. cs%id_pdevstress_2 > 0 .or. &
1205 cs%id_pstrainrate_1 > 0 .or. cs%id_pstrainrate_2 > 0) then
1206 call is_dynamics_post_data_2(cs, iss, g)
1207 endif
1208
1209 call disable_averaging(cs%diag)
1210end subroutine is_dynamics_post_data
1211
1212!> Calculate cell-centered, area-averaged, vertically integrated ice viscosity for diagnostics
1213subroutine ice_visc_diag(CS,G,ice_visc)
1214 type(ice_shelf_dyn_cs), intent(in) :: CS !< The ice shelf dynamics control structure
1215 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
1216 real, dimension(SZDI_(G),SZDJ_(G)), intent(out) :: ice_visc !< area-averaged vertically integrated ice viscosity
1217 !! [R L2 Z T-1 ~> Pa s m]
1218 integer :: i,j
1219
1220 ice_visc(:,:)=0.0
1221 if (cs%visc_qps==4) then
1222 do j=g%jsc,g%jec ; do i=g%isc,g%iec
1223 ice_visc(i,j) = (0.25 * g%IareaT(i,j)) * &
1224 ((cs%ice_visc(i,j,1) + cs%ice_visc(i,j,4)) + (cs%ice_visc(i,j,2) + cs%ice_visc(i,j,3)))
1225 enddo ; enddo
1226 else
1227 do j=g%jsc,g%jec ; do i=g%isc,g%iec
1228 ice_visc(i,j) = cs%ice_visc(i,j,1)*g%IareaT(i,j)
1229 enddo ; enddo
1230 endif
1231end subroutine ice_visc_diag
1232
1233!> Writes the total ice shelf kinetic energy and mass to an ascii file
1234subroutine write_ice_shelf_energy(CS, G, US, mass, area, day, time_step)
1235 type(ice_shelf_dyn_cs), intent(inout) :: cs !< The ice shelf dynamics control structure
1236 type(ocean_grid_type), intent(inout) :: g !< The grid structure used by the ice shelf.
1237 type(unit_scale_type), intent(in) :: us !< A structure containing unit conversion factors
1238 real, dimension(SZDI_(G),SZDJ_(G)), &
1239 intent(in) :: mass !< The mass per unit area of the ice shelf
1240 !! or sheet [R Z ~> kg m-2]
1241 real, dimension(SZDI_(G),SZDJ_(G)), &
1242 intent(in) :: area !< The ice shelf or ice sheet area [L2 ~> m2]
1243 type(time_type), intent(in) :: day !< The current model time.
1244 type(time_type), optional, intent(in) :: time_step !< The current time step
1245 ! Local variables
1246 type(time_type) :: dt ! A time_type version of the timestep.
1247 real, dimension(SZDI_(G),SZDJ_(G)) :: tmp1 ! A temporary array used in reproducing sums [various]
1248 real :: ke_tot ! The total kinetic energy [R Z L4 T-2 ~> J]
1249 real :: mass_tot ! The total mass [R Z L2 ~> kg]
1250 integer :: is, ie, js, je, isr, ier, jsr, jer, i, j
1251 character(len=32) :: mesg_intro, time_units, day_str, n_str, date_str
1252 integer :: start_of_day, num_days
1253 real :: reday ! Time in units given by CS%Timeunit, but often [days]
1254
1255 ! write_energy_time is the next integral multiple of energysavedays.
1256 if (present(time_step)) then
1257 dt = time_step
1258 else
1259 dt = set_time(seconds=2)
1260 endif
1261
1262 !CS%prev_IS_energy_calls tracks the ice sheet step, which is outputted in the energy file.
1263 if (cs%prev_IS_energy_calls == 0) then
1264 if (cs%energysave_geometric) then
1265 if (cs%energysavedays_geometric < cs%energysavedays) then
1266 cs%write_energy_time = day + cs%energysavedays_geometric
1267 cs%geometric_end_time = cs%Start_time + cs%energysavedays * &
1268 (1 + (day - cs%Start_time) / cs%energysavedays)
1269 else
1270 cs%write_energy_time = cs%Start_time + cs%energysavedays * &
1271 (1 + (day - cs%Start_time) / cs%energysavedays)
1272 endif
1273 else
1274 cs%write_energy_time = cs%Start_time + cs%energysavedays * &
1275 (1 + (day - cs%Start_time) / cs%energysavedays)
1276 endif
1277 elseif (day + (dt/2) <= cs%write_energy_time) then
1278 cs%prev_IS_energy_calls = cs%prev_IS_energy_calls + 1
1279 return ! Do not write this step
1280 else ! Determine the next write time before proceeding
1281 if (cs%energysave_geometric) then
1282 if (cs%write_energy_time + cs%energysavedays_geometric >= &
1283 cs%geometric_end_time) then
1284 cs%write_energy_time = cs%geometric_end_time
1285 cs%energysave_geometric = .false. ! stop geometric progression
1286 else
1287 cs%write_energy_time = cs%write_energy_time + cs%energysavedays_geometric
1288 endif
1289 cs%energysavedays_geometric = cs%energysavedays_geometric*2
1290 else
1291 cs%write_energy_time = cs%write_energy_time + cs%energysavedays
1292 endif
1293 endif
1294
1295 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
1296 isr = is - (g%isd-1) ; ier = ie - (g%isd-1) ; jsr = js - (g%jsd-1) ; jer = je - (g%jsd-1)
1297
1298 !calculate KE using cell-centered ice shelf velocity
1299 tmp1(:,:) = 0.0
1300 do j=js,je ; do i=is,ie
1301 tmp1(i,j) = 0.03125 * (mass(i,j) * area(i,j)) * &
1302 ((((cs%u_shelf(i-1,j-1)+cs%u_shelf(i,j))+(cs%u_shelf(i,j-1)+cs%u_shelf(i-1,j)))**2) + &
1303 (((cs%v_shelf(i-1,j-1)+cs%v_shelf(i,j))+(cs%v_shelf(i,j-1)+cs%v_shelf(i-1,j)))**2))
1304 enddo ; enddo
1305
1306 ke_tot = reproducing_sum(tmp1, isr, ier, jsr, jer, unscale=(us%RZL2_to_kg*us%L_T_to_m_s**2))
1307
1308 !calculate mass
1309 tmp1(:,:) = 0.0
1310 do j=js,je ; do i=is,ie
1311 tmp1(i,j) = mass(i,j) * area(i,j)
1312 enddo ; enddo
1313
1314 mass_tot = reproducing_sum(tmp1, isr, ier, jsr, jer, unscale=us%RZL2_to_kg)
1315
1316 if (is_root_pe()) then ! Only the root PE actually writes anything.
1317 if (day > cs%Start_time) then
1318 call open_ascii_file(cs%IS_fileenergy_ascii, trim(cs%IS_energyfile), action=append_file)
1319 else
1320 call open_ascii_file(cs%IS_fileenergy_ascii, trim(cs%IS_energyfile), action=writeonly_file)
1321 if (abs(cs%timeunit - 86400.0) < 1.0) then
1322 write(cs%IS_fileenergy_ascii,'(" Step,",7x,"Day,",8x,"Energy/Mass,",13x,"Total Mass")')
1323 write(cs%IS_fileenergy_ascii,'(12x,"[days]",10x,"[m2 s-2]",17x,"[kg]")')
1324 else
1325 if ((cs%timeunit >= 0.99) .and. (cs%timeunit < 1.01)) then
1326 time_units = " [seconds] "
1327 elseif ((cs%timeunit >= 3599.0) .and. (cs%timeunit < 3601.0)) then
1328 time_units = " [hours] "
1329 elseif ((cs%timeunit >= 86399.0) .and. (cs%timeunit < 86401.0)) then
1330 time_units = " [days] "
1331 elseif ((cs%timeunit >= 3.0e7) .and. (cs%timeunit < 3.2e7)) then
1332 time_units = " [years] "
1333 else
1334 write(time_units,'(9x,"[",es8.2," s] ")') cs%timeunit
1335 endif
1336
1337 write(cs%IS_fileenergy_ascii,'(" Step,",7x,"Time,",7x,"Energy/Mass,",13x,"Total Mass")')
1338 write(cs%IS_fileenergy_ascii,'(A25,3x,"[m2 s-2]",17x,"[kg]")') time_units
1339 endif
1340 endif
1341
1342 call get_time(day, start_of_day, num_days)
1343
1344 if (abs(cs%timeunit - 86400.0) < 1.0) then
1345 reday = real(num_days)+ (real(start_of_day)/86400.0)
1346 else
1347 reday = real(num_days)*(86400.0/cs%timeunit) + real(start_of_day)/abs(cs%timeunit)
1348 endif
1349
1350 if (reday < 1.0e8) then ; write(day_str, '(F12.3)') reday
1351 elseif (reday < 1.0e11) then ; write(day_str, '(F15.3)') reday
1352 else ; write(day_str, '(ES15.9)') reday ; endif
1353
1354 if (cs%prev_IS_energy_calls < 1000000) then ; write(n_str, '(I6)') cs%prev_IS_energy_calls
1355 else ; write(n_str, '(I0)') cs%prev_IS_energy_calls ; endif
1356
1357 write(cs%IS_fileenergy_ascii,'(A,",",A,", En ",ES22.16,", M ",ES11.5)') &
1358 trim(n_str), trim(day_str), us%L_T_to_m_s**2*ke_tot/mass_tot, us%RZL2_to_kg*mass_tot
1359 endif
1360
1361 cs%prev_IS_energy_calls = cs%prev_IS_energy_calls + 1
1362end subroutine write_ice_shelf_energy
1363
1364!> This subroutine takes the velocity (on the Bgrid) and timesteps h_t = - div (uh) once.
1365!! Additionally, it will update the volume of ice in partially-filled cells, and update
1366!! hmask accordingly
1367subroutine ice_shelf_advect(CS, ISS, G, time_step, Time, calve_ice_shelf_bergs)
1368 type(ice_shelf_dyn_cs), intent(inout) :: CS !< The ice shelf dynamics control structure
1369 type(ice_shelf_state), intent(inout) :: ISS !< A structure with elements that describe
1370 !! the ice-shelf state
1371 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
1372 real, intent(in) :: time_step !< time step [T ~> s]
1373 type(time_type), intent(in) :: Time !< The current model time
1374 logical, intent(in) :: calve_ice_shelf_bergs !< If true, track ice shelf flux through a
1375 !! static ice shelf, so that it can be converted into icebergs
1376
1377! 3/8/11 DNG
1378!
1379! This subroutine takes the velocity (on the Bgrid) and timesteps h_t = - div (uh) once.
1380! ADDITIONALLY, it will update the volume of ice in partially-filled cells, and update
1381! hmask accordingly
1382!
1383! The flux overflows are included here. That is because they will be used to advect 3D scalars
1384! into partial cells
1385
1386 real, dimension(SZDI_(G),SZDJ_(G)) :: h_after_flux1, h_after_flux2 ! Ice thicknesses [Z ~> m].
1387 real, dimension(SZDIB_(G),SZDJ_(G)) :: uh_ice ! The accumulated zonal ice volume flux [Z L2 ~> m3]
1388 real, dimension(SZDI_(G),SZDJB_(G)) :: vh_ice ! The accumulated meridional ice volume flux [Z L2 ~> m3]
1389 type(loop_bounds_type) :: LB
1390 integer :: isd, ied, jsd, jed, i, j, isc, iec, jsc, jec, stencil
1391
1392 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1393 isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
1394
1395 uh_ice(:,:) = 0.0
1396 vh_ice(:,:) = 0.0
1397
1398 h_after_flux1(:,:) = 0.0
1399 h_after_flux2(:,:) = 0.0
1400 ! call MOM_mesg("MOM_ice_shelf.F90: ice_shelf_advect called")
1401
1402 do j=jsd,jed ; do i=isd,ied ; if (cs%h_bdry_val(i,j) /= 0.0) then
1403 iss%h_shelf(i,j) = cs%h_bdry_val(i,j)
1404 endif ; enddo ; enddo
1405
1406 stencil = 2
1407 if (modulo(cs%first_direction_IS,2)==0) then
1408 !x first
1409 lb%ish = g%isc ; lb%ieh = g%iec ; lb%jsh = g%jsc-stencil ; lb%jeh = g%jec+stencil
1410 if (lb%jsh < jsd) call mom_error(fatal, &
1411 "ice_shelf_advect: Halo is too small for the ice thickness advection stencil.")
1412 call ice_shelf_advect_thickness_x(cs, g, lb, time_step, iss%hmask, iss%h_shelf, h_after_flux1, uh_ice)
1413 call pass_var(h_after_flux1, g%domain)
1414 lb%ish = g%isc ; lb%ieh = g%iec ; lb%jsh = g%jsc ; lb%jeh = g%jec
1415 call ice_shelf_advect_thickness_y(cs, g, lb, time_step, iss%hmask, h_after_flux1, h_after_flux2, vh_ice)
1416 else
1417 ! y first
1418 lb%ish = g%isc-stencil ; lb%ieh = g%iec+stencil ; lb%jsh = g%jsc ; lb%jeh = g%jec
1419 if (lb%ish < isd) call mom_error(fatal, &
1420 "ice_shelf_advect: Halo is too small for the ice thickness advection stencil.")
1421 call ice_shelf_advect_thickness_y(cs, g, lb, time_step, iss%hmask, iss%h_shelf, h_after_flux1, vh_ice)
1422 call pass_var(h_after_flux1, g%domain)
1423 lb%ish = g%isc ; lb%ieh = g%iec ; lb%jsh = g%jsc ; lb%jeh = g%jec
1424 call ice_shelf_advect_thickness_x(cs, g, lb, time_step, iss%hmask, h_after_flux1, h_after_flux2, uh_ice)
1425 endif
1426 call pass_var(h_after_flux2, g%domain)
1427
1428 do j=jsd,jed
1429 do i=isd,ied
1430 if (iss%hmask(i,j) == 1) iss%h_shelf(i,j) = h_after_flux2(i,j)
1431 enddo
1432 enddo
1433
1434 if (cs%moving_shelf_front) then
1435 call shelf_advance_front(cs, iss, g, iss%hmask, uh_ice, vh_ice)
1436 if (cs%min_thickness_simple_calve > 0.0) then
1437 call ice_shelf_min_thickness_calve(g, iss%h_shelf, iss%area_shelf_h, iss%hmask, &
1438 cs%min_thickness_simple_calve)
1439 endif
1440 if (cs%calve_to_mask) then
1441 call calve_to_mask(g, iss%h_shelf, iss%area_shelf_h, iss%hmask, cs%calve_mask)
1442 endif
1443 elseif (calve_ice_shelf_bergs) then
1444 !advect the front to create partially-filled cells
1445 call shelf_advance_front(cs, iss, g, iss%hmask, uh_ice, vh_ice)
1446 !add mass of the partially-filled cells to calving field, which is used to initialize icebergs
1447 !Then, remove the partially-filled cells from the ice shelf
1448 iss%calving(:,:) = 0.0
1449 iss%calving_hflx(:,:) = 0.0
1450 do j=jsc,jec ; do i=isc,iec
1451 if (iss%hmask(i,j)==2) then
1452 iss%calving(i,j) = (iss%h_shelf(i,j) * cs%density_ice) * &
1453 (iss%area_shelf_h(i,j) * g%IareaT(i,j)) / time_step
1454 iss%calving_hflx(i,j) = (cs%Cp_ice * cs%t_shelf(i,j)) * &
1455 ((iss%h_shelf(i,j) * cs%density_ice) * &
1456 (iss%area_shelf_h(i,j) * g%IareaT(i,j)))
1457 iss%h_shelf(i,j) = 0.0 ; iss%area_shelf_h(i,j) = 0.0 ; iss%hmask(i,j) = 0.0
1458 endif
1459 enddo ; enddo
1460 endif
1461
1462 do j=jsc,jec ; do i=isc,iec
1463 iss%mass_shelf(i,j) = iss%h_shelf(i,j) * cs%density_ice
1464 enddo ; enddo
1465
1466 call pass_var(iss%mass_shelf, g%domain, complete=.false.)
1467 call pass_var(iss%h_shelf, g%domain, complete=.false.)
1468 call pass_var(iss%area_shelf_h, g%domain, complete=.false.)
1469 call pass_var(iss%hmask, g%domain, complete=.true.)
1470
1471 call update_velocity_masks(cs, g, iss%hmask, cs%umask, cs%vmask, cs%u_face_mask, cs%v_face_mask)
1472
1473end subroutine ice_shelf_advect
1474
1475!>This subroutine computes u- and v-velocities of the ice shelf iterating on non-linear ice viscosity
1476!subroutine ice_shelf_solve_outer(CS, ISS, G, US, u_shlf, v_shlf, iters, time)
1477subroutine ice_shelf_solve_outer(CS, ISS, G, US, u_shlf, v_shlf, taudx, taudy, iters, Time)
1478 type(ice_shelf_dyn_cs), intent(inout) :: CS !< The ice shelf dynamics control structure
1479 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
1480 !! the ice-shelf state
1481 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
1482 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
1483 real, dimension(SZDIB_(G),SZDJB_(G)), &
1484 intent(inout) :: u_shlf !< The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
1485 real, dimension(SZDIB_(G),SZDJB_(G)), &
1486 intent(inout) :: v_shlf !< The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
1487 integer, intent(out) :: iters !< The number of iterations used in the solver.
1488 type(time_type), intent(in) :: Time !< The current model time
1489
1490 real, dimension(SZDIB_(G),SZDJB_(G)), &
1491 intent(out) :: taudx !< Driving x-stress at q-points [R L3 Z T-2 ~> kg m s-2]
1492 real, dimension(SZDIB_(G),SZDJB_(G)), &
1493 intent(out) :: taudy !< Driving y-stress at q-points [R L3 Z T-2 ~> kg m s-2]
1494 !real, dimension(SZDIB_(G),SZDJB_(G)) :: u_bdry_cont ! Boundary u-stress contribution [R L3 Z T-2 ~> kg m s-2]
1495 !real, dimension(SZDIB_(G),SZDJB_(G)) :: v_bdry_cont ! Boundary v-stress contribution [R L3 Z T-2 ~> kg m s-2]
1496 real, dimension(SZDIB_(G),SZDJB_(G)) :: Au, Av ! The retarding lateral stress contributions [R L3 Z T-2 ~> kg m s-2]
1497 real, dimension(SZDIB_(G),SZDJB_(G)) :: u_last, v_last ! Previous velocities [L T-1 ~> m s-1]
1498 real, dimension(SZDIB_(G),SZDJB_(G)) :: H_node ! Ice shelf thickness at corners [Z ~> m].
1499 real, dimension(SZDI_(G),SZDJ_(G)) :: float_cond ! If GL_regularize=true, indicates cells containing
1500 ! the grounding line (float_cond=1) or not (float_cond=0)
1501 real, dimension(SZDIB_(G),SZDJB_(G)) :: Normvec ! Velocities used for convergence [L2 T-2 ~> m2 s-2]
1502 character(len=160) :: mesg ! The text of an error message
1503 integer :: conv_flag, i, j, k,l, iter, nodefloat
1504 integer :: Isdq, Iedq, Jsdq, Jedq, isd, ied, jsd, jed
1505 integer :: Iscq, Iecq, Jscq, Jecq, isc, iec, jsc, jec
1506 real :: err_max, err_tempu, err_tempv, err_init ! Errors in [R L3 Z T-2 ~> kg m s-2] or [L T-1 ~> m s-1]
1507 real :: ew_prev_err ! Previous outer residual for Eisenstat-Walker CG tolerance (same units as err_max)
1508 real :: max_vel ! The maximum velocity magnitude [L T-1 ~> m s-1]
1509 real :: tempu, tempv ! Temporary variables with velocity magnitudes [L T-1 ~> m s-1]
1510 real :: Norm, PrevNorm ! Velocities used to assess convergence [L T-1 ~> m s-1]
1511 real :: rhoi_rhow ! The density of ice divided by a typical water density [nondim]
1512 integer :: Is_sum, Js_sum, Ie_sum, Je_sum ! Loop bounds for global sums or arrays starting at 1.
1513 integer :: Iscq_sv, Jscq_sv ! Starting loop bound for sum_vec
1514
1515 isdq = g%IsdB ; iedq = g%IedB ; jsdq = g%JsdB ; jedq = g%JedB
1516 iscq = g%IscB ; iecq = g%IecB ; jscq = g%JscB ; jecq = g%JecB
1517 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1518 isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
1519 rhoi_rhow = cs%density_ice / cs%density_ocean_avg
1520
1521 taudx(:,:) = 0.0 ; taudy(:,:) = 0.0
1522 au(:,:) = 0.0 ; av(:,:) = 0.0
1523
1524 ! need to make these conditional on GL interpolation
1525 cs%float_cond(:,:) = 0.0 ; h_node(:,:) = 0.0
1526 !CS%ground_frac(:,:) = 0.0
1527
1528 if (.not. cs%GL_couple) then
1529 do j=g%jsc,g%jec ; do i=g%isc,g%iec
1530 if (rhoi_rhow * max(iss%h_shelf(i,j),cs%min_h_shelf) - cs%bed_elev(i,j) > 0) then
1531 cs%ground_frac(i,j) = 1.0
1532 cs%OD_av(i,j) =0.0
1533 endif
1534 enddo ; enddo
1535 endif
1536
1537 call calc_shelf_driving_stress(cs, iss, g, us, taudx, taudy, cs%OD_av)
1538 call pass_vector(taudx, taudy, g%domain, to_all, bgrid_ne)
1539 ! This is to determine which cells contain the grounding line, the criterion being that the cell
1540 ! is ice-covered, with some nodes floating and some grounded flotation condition is estimated by
1541 ! assuming topography is cellwise constant and H is bilinear in a cell; floating where
1542 ! rho_i/rho_w * H_node - D is negative
1543
1544 ! need to make this conditional on GL interp
1545
1546 if (cs%GL_regularize) then
1547
1548 call interpolate_h_to_b(g, iss%h_shelf, iss%hmask, h_node, cs%min_h_shelf)
1549
1550 do j=g%jsc,g%jec ; do i=g%isc,g%iec
1551 nodefloat = 0
1552
1553 do l=0,1 ; do k=0,1
1554 if ((iss%hmask(i,j) == 1 .or. iss%hmask(i,j)==3) .and. &
1555 (rhoi_rhow * h_node(i-1+k,j-1+l) - cs%bed_elev(i,j) <= 0)) then
1556 nodefloat = nodefloat + 1
1557 endif
1558 enddo ; enddo
1559 if ((nodefloat > 0) .and. (nodefloat < 4)) then
1560 cs%float_cond(i,j) = 1.0
1561 cs%ground_frac(i,j) = 1.0
1562 endif
1563 enddo ; enddo
1564
1565 call pass_var(cs%float_cond, g%Domain, complete=.false.)
1566 call pass_var(cs%ground_frac, g%domain, complete=.false.)
1567
1568 endif
1569
1570 call calc_shelf_taub(cs, iss, g, us, u_shlf, v_shlf)
1571 call pass_var(cs%basal_traction, g%domain, complete=.true.)
1572 call calc_shelf_visc(cs, iss, g, us, u_shlf, v_shlf)
1573 call pass_var(cs%ice_visc, g%domain)
1574
1575 ! This makes sure basal stress is only applied when it is supposed to be
1576 if (cs%GL_regularize) then
1577 do j=g%jsd,g%jed ; do i=g%isd,g%ied
1578 if (cs%ground_frac(i,j)/=1.0) cs%basal_traction(i,j) = 0.0
1579 enddo ; enddo
1580 else
1581 do j=g%jsd,g%jed ; do i=g%isd,g%ied
1582 cs%basal_traction(i,j) = cs%basal_traction(i,j) * cs%ground_frac(i,j)
1583 enddo ; enddo
1584 endif
1585
1586 if (cs%nonlin_solve_err_mode == 1) then
1587
1588 au(:,:) = 0.0 ; av(:,:) = 0.0
1589
1590 call cg_action(cs, au, av, u_shlf, v_shlf, cs%Phi, cs%Phisub, cs%umask, cs%vmask, iss%hmask, h_node, &
1591 cs%ice_visc, cs%float_cond, cs%bed_elev, cs%basal_traction, &
1592 g, us, g%isc-1, g%iec+1, g%jsc-1, g%jec+1, rhoi_rhow, use_newton_in=.false.)
1593 call pass_vector(au, av, g%domain, to_all, bgrid_ne)
1594
1595 err_init = 0 ; err_tempu = 0 ; err_tempv = 0
1596 do j=g%JscB,g%JecB ; do i=g%IscB,g%IecB
1597 if (cs%umask(i,j) == 1) then
1598 err_tempu = abs(au(i,j) - taudx(i,j))
1599 if (err_tempu >= err_init) err_init = err_tempu
1600 endif
1601 if (cs%vmask(i,j) == 1) then
1602 err_tempv = abs(av(i,j) - taudy(i,j))
1603 if (err_tempv >= err_init) err_init = err_tempv
1604 endif
1605 enddo ; enddo
1606
1607 call max_across_pes(err_init)
1608 elseif (cs%nonlin_solve_err_mode == 3) then
1609 normvec(:,:) = 0.0
1610
1611 ! Determine the loop limits for sums, bearing in mind that the arrays will be starting at 1.
1612 ! Includes the edge of the tile is at the western/southern bdry (if symmetric)
1613 if ((isc+g%idg_offset==g%isg) .and. (.not. cs%reentrant_x)) then
1614 is_sum = iscq + (1-isdq) ; iscq_sv = iscq
1615 else
1616 is_sum = isc + (1-isdq) ; iscq_sv = isc
1617 endif
1618 if ((jsc+g%jdg_offset==g%jsg) .and. (.not. cs%reentrant_y)) then
1619 js_sum = jscq + (1-jsdq) ; jscq_sv = jscq
1620 else
1621 js_sum = jsc + (1-jsdq) ; jscq_sv = jsc
1622 endif
1623 ie_sum = iecq + (1-isdq) ; je_sum = jecq + (1-jsdq)
1624
1625 do j=jscq_sv,jecq ; do i=iscq_sv,iecq
1626 if (cs%umask(i,j) == 1) normvec(i,j) = normvec(i,j) + u_shlf(i,j)**2
1627 if (cs%vmask(i,j) == 1) normvec(i,j) = normvec(i,j) + v_shlf(i,j)**2
1628 enddo ; enddo
1629 norm = sqrt( reproducing_sum( normvec, is_sum, ie_sum, js_sum, je_sum, unscale=us%L_T_to_m_s**2 ) )
1630 endif
1631
1632 u_last(:,:) = u_shlf(:,:) ; v_last(:,:) = v_shlf(:,:)
1633 cs%cg_tol_newton = cs%cg_tolerance
1634 ew_prev_err = err_init
1635
1636 !! begin loop
1637
1638 do iter=1,50
1639
1640 call ice_shelf_solve_inner(cs, iss, g, us, u_shlf, v_shlf, taudx, taudy, h_node, cs%float_cond, &
1641 iss%hmask, conv_flag, iters, time, cs%Phi, cs%Phisub)
1642
1643 if (cs%debug) then
1644 call qchksum(u_shlf, "u shelf", g%HI, haloshift=2, unscale=us%L_T_to_m_s)
1645 call qchksum(v_shlf, "v shelf", g%HI, haloshift=2, unscale=us%L_T_to_m_s)
1646 endif
1647
1648 write(mesg,*) "ice_shelf_solve_outer: linear solve done in ",iters," iterations"
1649 call mom_mesg(mesg, 5)
1650
1651 call calc_shelf_taub(cs, iss, g, us, u_shlf, v_shlf)
1652 call pass_var(cs%basal_traction, g%domain, complete=.true.)
1653 call calc_shelf_visc(cs, iss, g, us, u_shlf, v_shlf)
1654 call pass_var(cs%ice_visc, g%domain, complete=.false.)
1655 call pass_var(cs%newton_str_sh, g%domain, complete=.false.)
1656 call pass_var(cs%newton_visc_factor, g%domain, complete=.true.)
1657 call pass_var(cs%newton_drag_coef, g%domain)
1658 call pass_vector(cs%newton_str_ux, cs%newton_str_vy, g%domain, to_all, agrid)
1659 call pass_vector(cs%newton_umid, cs%newton_vmid, g%domain, to_all, agrid)
1660
1661 ! makes sure basal stress is only applied when it is supposed to be
1662 if (cs%GL_regularize) then
1663 do j=g%jsd,g%jed ; do i=g%isd,g%ied
1664 if (cs%ground_frac(i,j)/=1.0) cs%basal_traction(i,j) = 0.0
1665 enddo ; enddo
1666 else
1667 do j=g%jsd,g%jed ; do i=g%isd,g%ied
1668 cs%basal_traction(i,j) = cs%basal_traction(i,j) * cs%ground_frac(i,j)
1669 enddo ; enddo
1670 endif
1671
1672 if (cs%nonlin_solve_err_mode == 1) then
1673
1674 au(:,:) = 0 ; av(:,:) = 0
1675
1676 call cg_action(cs, au, av, u_shlf, v_shlf, cs%Phi, cs%Phisub, cs%umask, cs%vmask, iss%hmask, h_node, &
1677 cs%ice_visc, cs%float_cond, cs%bed_elev, cs%basal_traction, &
1678 g, us, g%isc-1, g%iec+1, g%jsc-1, g%jec+1, rhoi_rhow, use_newton_in=.false.)
1679
1680 call pass_vector(au, av, g%domain, to_all, bgrid_ne)
1681
1682 err_max = 0
1683
1684 do j=g%jscB,g%jecB ; do i=g%iscB,g%iecB
1685 if (cs%umask(i,j) == 1) then
1686 err_tempu = abs(au(i,j) - taudx(i,j))
1687 if (err_tempu >= err_max) err_max = err_tempu
1688 endif
1689 if (cs%vmask(i,j) == 1) then
1690 err_tempv = abs(av(i,j) - taudy(i,j))
1691 if (err_tempv >= err_max) err_max = err_tempv
1692 endif
1693 enddo ; enddo
1694
1695 call max_across_pes(err_max)
1696
1697 elseif (cs%nonlin_solve_err_mode == 2) then
1698
1699 err_max=0. ; max_vel = 0 ; tempu = 0 ; tempv = 0 ; err_tempu = 0
1700 do j=g%jscB,g%jecB ; do i=g%iscB,g%iecB
1701 if (cs%umask(i,j) == 1) then
1702 err_tempu = abs(u_last(i,j)-u_shlf(i,j))
1703 if (err_tempu >= err_max) err_max = err_tempu
1704 tempu = u_shlf(i,j)
1705 else
1706 tempu = 0.0
1707 endif
1708 if (cs%vmask(i,j) == 1) then
1709 err_tempv = max(abs(v_last(i,j)-v_shlf(i,j)), err_tempu)
1710 if (err_tempv >= err_max) err_max = err_tempv
1711 tempv = sqrt((v_shlf(i,j)**2) + (tempu**2))
1712 endif
1713 if (tempv >= max_vel) max_vel = tempv
1714 enddo ; enddo
1715
1716 u_last(:,:) = u_shlf(:,:)
1717 v_last(:,:) = v_shlf(:,:)
1718
1719 call max_across_pes(max_vel)
1720 call max_across_pes(err_max)
1721 err_init = max_vel
1722
1723 elseif (cs%nonlin_solve_err_mode == 3) then
1724 prevnorm = norm ; norm = 0.0 ; normvec=0.0
1725 do j=jscq_sv,jecq ; do i=iscq_sv,iecq
1726 if (cs%umask(i,j) == 1) normvec(i,j) = normvec(i,j) + u_shlf(i,j)**2
1727 if (cs%vmask(i,j) == 1) normvec(i,j) = normvec(i,j) + v_shlf(i,j)**2
1728 enddo ; enddo
1729 norm = sqrt( reproducing_sum( normvec, is_sum, ie_sum, js_sum, je_sum, unscale=us%L_T_to_m_s**2 ) )
1730 err_max = 2.*abs(norm-prevnorm) ; err_init = norm+prevnorm
1731 endif
1732
1733 write(mesg,*) "ice_shelf_solve_outer: nonlinear fractional residual = ", err_max/err_init
1734 call mom_mesg(mesg, 5)
1735
1736 if (err_max <= cs%newton_after_tolerance * err_init .and. .not. cs%doing_newton) then
1737 cs%doing_newton = .true.
1738 ew_prev_err = err_max ! seed Eisenstat-Walker with residual at the Newton switch point
1739 write(mesg,*) "ice_shelf_solve_outer: switching to Newton iterations at iter = ", iter
1740 call mom_mesg(mesg, 5)
1741 endif
1742
1743 ! Eisenstat-Walker Choice II (Eisenstat & Walker 1994): η_k = γ*(||F_k||/||F_{k-1}||)^α
1744 ! with γ=0.9, α=2. Uses the ratio of consecutive outer residuals so that the CG
1745 ! tolerance scales linearly with the current error (enabling quadratic outer convergence)
1746 ! without over-tightening at later Newton steps. The first Newton step uses the standard
1747 ! cg_tolerance (ratio = 1 on entry).
1748 if (cs%doing_newton .and. cs%newton_adapt_cg_tol) then
1749 cs%cg_tol_newton = min(cs%cg_tolerance, 0.9 * (err_max / ew_prev_err)**2)
1750 ew_prev_err = err_max
1751 endif
1752
1753 if (err_max <= cs%nonlinear_tolerance * err_init) then
1754 exit
1755 endif
1756
1757 enddo
1758 cs%doing_newton = .false.
1759 cs%cg_tol_newton = cs%cg_tolerance
1760
1761 write(mesg,*) "ice_shelf_solve_outer: nonlinear fractional residual = ", err_max/err_init
1762 call mom_mesg(mesg)
1763 write(mesg,*) "ice_shelf_solve_outer: exiting nonlinear solve after ",iter," iterations"
1764 call mom_mesg(mesg)
1765
1766end subroutine ice_shelf_solve_outer
1767
1768subroutine ice_shelf_solve_inner(CS, ISS, G, US, u_shlf, v_shlf, taudx, taudy, H_node, float_cond, &
1769 hmask, conv_flag, iters, time, Phi, Phisub)
1770 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
1771 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
1772 !! the ice-shelf state
1773 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
1774 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
1775 real, dimension(SZDIB_(G),SZDJB_(G)), &
1776 intent(inout) :: u_shlf !< The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
1777 real, dimension(SZDIB_(G),SZDJB_(G)), &
1778 intent(inout) :: v_shlf !< The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
1779 real, dimension(SZDIB_(G),SZDJB_(G)), &
1780 intent(in) :: taudx !< The x-direction driving stress [R L3 Z T-2 ~> kg m s-2]
1781 real, dimension(SZDIB_(G),SZDJB_(G)), &
1782 intent(in) :: taudy !< The y-direction driving stress [R L3 Z T-2 ~> kg m s-2]
1783 real, dimension(SZDIB_(G),SZDJB_(G)), &
1784 intent(in) :: H_node !< The ice shelf thickness at nodal (corner)
1785 !! points [Z ~> m].
1786 real, dimension(SZDI_(G),SZDJ_(G)), &
1787 intent(in) :: float_cond !< If GL_regularize=true, indicates cells containing
1788 !! the grounding line (float_cond=1) or not (float_cond=0)
1789 real, dimension(SZDI_(G),SZDJ_(G)), &
1790 intent(in) :: hmask !< A mask indicating which tracer points are
1791 !! partly or fully covered by an ice-shelf
1792 integer, intent(out) :: conv_flag !< A flag indicating whether (1) or not (0) the
1793 !! iterations have converged to the specified tolerance
1794 integer, intent(out) :: iters !< The number of iterations used in the solver.
1795 type(time_type), intent(in) :: Time !< The current model time
1796 real, dimension(8,4,SZDI_(G),SZDJ_(G)), &
1797 intent(in) :: Phi !< The gradients of bilinear basis elements at Gaussian
1798 !! quadrature points surrounding the cell vertices [L-1 ~> m-1].
1799 real, dimension(:,:,:,:,:,:), &
1800 intent(in) :: Phisub !< Quadrature structure weights at subgridscale
1801 !! locations for finite element calculations [nondim]
1802! one linear solve (nonlinear iteration) of the solution for velocity
1803
1804! in this subroutine:
1805! RHS = taud
1806! diagonal of matrix is found (for Jacobi precondition)
1807! CG iteration is carried out for max. iterations or until convergence
1808
1809! assumed - u, v, taud, visc, basal_traction are valid on the halo
1810
1811 real, dimension(SZDIB_(G),SZDJB_(G)) :: &
1812 Ru, Rv, & ! Residuals in the stress calculations [R L3 Z T-2 ~> m kg s-2]
1813 Ru_old, Rv_old, & ! Previous values of Ru and Rv [R L3 Z T-2 ~> m kg s-2]
1814 Zu, Zv, & ! Contributions to velocity changes [L T-1 ~> m s-1]
1815 Zu_old, Zv_old, & ! Previous values of Zu and Zv [L T-1 ~> m s-1]
1816 DIAGu, DIAGv, & ! Diagonals with units like Ru/Zu [R L2 Z T-1 ~> kg s-1]
1817 RHSu, RHSv, & ! Right hand side of the stress balance [R L3 Z T-2 ~> m kg s-2]
1818 Au, Av, & ! The retarding lateral stress contributions [R L3 Z T-2 ~> kg m s-2]
1819 Du, Dv, & ! Velocity changes [L T-1 ~> m s-1]
1820 sum_vec ! Sum of squares of residuals in stress calculations [m2 kg2 s-4]
1821 real, dimension(SZDIB_(G),SZDJB_(G),3) :: sum_vec_3d ! Array used for various residuals
1822 ! sum_vec_3d(:,:,1:2) [m s-1] [m kg s-2]
1823 ! sum_vec_3d(:,:,3) [m2 kg2 s-4]
1824 real :: beta_k ! Ratio of residuals used to update search direction [nondim]
1825 real :: resid0tol2 ! Convergence tolerance times the initial residual [m2 kg2 s-4]
1826 real :: sv3dsum ! An unused variable returned when taking global sum of residuals [various]
1827 real :: sv3dsums(3) ! The index-wise global sums of sum_vec_3d
1828 ! sv3dsums(:,:,1:2) [m s-1] [m kg s-2]
1829 ! sv3dsums(:,:,3) [m2 kg2 s-4]
1830 real :: alpha_k ! A scaling factor for iterative corrections [nondim]
1831 real :: resid_scale ! A scaling factor for redimensionalizing the global residuals
1832 ! [L T-1 ~> m s-1] [R L3 Z T-2 ~> m kg s-2]
1833 real :: resid2_scale ! A scaling factor for redimensionalizing the global squared residuals
1834 ! [R2 L6 Z2 T-4 ~> m2 kg2 s-4]
1835 real :: rhoi_rhow ! The density of ice divided by a typical water density [nondim]
1836 integer :: cg_halo ! Number of halo vertices to include during a CG iteration
1837 integer :: max_cg_halo ! Maximum possible number of halo vertices to include in the CG iterations
1838 integer :: iter, i, j, isd, ied, jsd, jed, isc, iec, jsc, jec, is, js, ie, je
1839 integer :: Is_sum, Js_sum, Ie_sum, Je_sum ! Loop bounds for global sums or arrays starting at 1.
1840 integer :: Isdq, Iedq, Jsdq, Jedq, Iscq, Iecq, Jscq, Jecq, nx_halo, ny_halo
1841 integer :: Iscq_sv, Jscq_sv ! Starting loop bound for sum_vec
1842
1843 isdq = g%IsdB ; iedq = g%IedB ; jsdq = g%JsdB ; jedq = g%JedB
1844 iscq = g%IscB ; iecq = g%IecB ; jscq = g%JscB ; jecq = g%JecB
1845 ny_halo = g%domain%njhalo ; nx_halo = g%domain%nihalo
1846 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
1847 isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
1848
1849 rhoi_rhow = cs%density_ice / cs%density_ocean_avg
1850
1851 zu(:,:) = 0 ; zv(:,:) = 0 ; diagu(:,:) = 0 ; diagv(:,:) = 0
1852 ru(:,:) = 0 ; rv(:,:) = 0 ; au(:,:) = 0 ; av(:,:) = 0 ; rhsu(:,:) = 0 ; rhsv(:,:) = 0
1853 du(:,:) = 0 ; dv(:,:) = 0
1854
1855 ! Determine the loop limits for sums, bearing in mind that the arrays will be starting at 1.
1856 ! Includes the edge of the tile is at the western/southern bdry (if symmetric)
1857 if ((isc+g%idg_offset==g%isg) .and. (.not. cs%reentrant_x)) then
1858 is_sum = iscq + (1-isdq) ; iscq_sv = iscq
1859 else
1860 is_sum = isc + (1-isdq) ; iscq_sv = isc
1861 endif
1862 if ((jsc+g%jdg_offset==g%jsg) .and. (.not. cs%reentrant_y)) then
1863 js_sum = jscq + (1-jsdq) ; jscq_sv = jscq
1864 else
1865 js_sum = jsc + (1-jsdq) ; jscq_sv = jsc
1866 endif
1867 ie_sum = iecq + (1-isdq) ; je_sum = jecq + (1-jsdq)
1868
1869 rhsu(:,:) = taudx(:,:) ; rhsv(:,:) = taudy(:,:)
1870
1871 call pass_vector(rhsu, rhsv, g%domain, to_all, bgrid_ne, complete=.false.)
1872
1873 call matrix_diagonal(cs, g, us, float_cond, h_node, cs%ice_visc, cs%basal_traction, &
1874 hmask, rhoi_rhow, phi, phisub, diagu, diagv)
1875
1876 call pass_vector(diagu, diagv, g%domain, to_all, bgrid_ne, complete=.false.)
1877
1878 call cg_action(cs, au, av, u_shlf, v_shlf, phi, phisub, cs%umask, cs%vmask, hmask, &
1879 h_node, cs%ice_visc, float_cond, cs%bed_elev, cs%basal_traction, &
1880 g, us, isc-1, iec+1, jsc-1, jec+1, rhoi_rhow, use_newton_in=.false.)
1881
1882 call pass_vector(au, av, g%domain, to_all, bgrid_ne, complete=.true.)
1883
1884 ru(:,:) = (rhsu(:,:) - au(:,:)) ; rv(:,:) = (rhsv(:,:) - av(:,:))
1885 resid_scale = us%s_to_T*(us%RZL2_to_kg*us%L_T_to_m_s**2)
1886 resid2_scale = ((us%RZ_to_kg_m2*us%L_to_m)*us%L_T_to_m_s**2)**2
1887
1888 sum_vec(:,:) = 0.0
1889 do j=jscq_sv,jecq ; do i=iscq_sv,iecq
1890 if (cs%umask(i,j) == 1) sum_vec(i,j) = resid2_scale*ru(i,j)**2
1891 if (cs%vmask(i,j) == 1) sum_vec(i,j) = sum_vec(i,j) + resid2_scale*rv(i,j)**2
1892 enddo ; enddo
1893
1894 !resid0 = sqrt(reproducing_sum( sum_vec, Is_sum, Ie_sum, Js_sum, Je_sum ))
1895 resid0tol2 = cs%cg_tol_newton**2 * reproducing_sum( sum_vec, is_sum, ie_sum, js_sum, je_sum )
1896
1897 do j=jsdq,jedq ; do i=isdq,iedq
1898 if (cs%umask(i,j) == 1 .AND.(diagu(i,j)/=0)) zu(i,j) = ru(i,j) / diagu(i,j)
1899 if (cs%vmask(i,j) == 1 .AND.(diagv(i,j)/=0)) zv(i,j) = rv(i,j) / diagv(i,j)
1900 enddo ; enddo
1901
1902 du(:,:) = zu(:,:) ; dv(:,:) = zv(:,:)
1903
1904 if (g%symmetric) then
1905 max_cg_halo=min(nx_halo,ny_halo)
1906 else
1907 max_cg_halo=min(nx_halo,ny_halo)-1
1908 endif
1909 cg_halo = max_cg_halo
1910 conv_flag = 0
1911
1912 !!!!!!!!!!!!!!!!!!
1913 !! !!
1914 !! MAIN CG LOOP !!
1915 !! !!
1916 !!!!!!!!!!!!!!!!!!
1917
1918 ! initially, c-grid data is valid up to 3 halo nodes out
1919
1920 do iter = 1,cs%cg_max_iterations
1921
1922 ! we can never assume that any arrays are legit more than 3 vertices past
1923 ! the computational domain - this is their state in the initial iteration
1924
1925 is = isc - cg_halo ; ie = iecq + cg_halo
1926 js = jsc - cg_halo ; je = jecq + cg_halo
1927
1928 au(:,:) = 0 ; av(:,:) = 0
1929
1930 call cg_action(cs, au, av, du, dv, phi, phisub, cs%umask, cs%vmask, hmask, &
1931 h_node, cs%ice_visc, float_cond, cs%bed_elev, cs%basal_traction, &
1932 g, us, is, ie, js, je, rhoi_rhow)
1933
1934 ! Au, Av valid region moves in by 1
1935
1936 call pass_vector(au,av,g%domain, to_all, bgrid_ne)
1937
1938 sum_vec_3d(:,:,1:2) = 0.0
1939
1940 do j=jscq_sv,jecq ; do i=iscq_sv,iecq
1941 if (cs%umask(i,j) == 1) then
1942 sum_vec_3d(i,j,1) = resid_scale * (zu(i,j) * ru(i,j))
1943 sum_vec_3d(i,j,2) = resid_scale * (du(i,j) * au(i,j))
1944 ru_old(i,j) = ru(i,j) ; zu_old(i,j) = zu(i,j)
1945 endif
1946 if (cs%vmask(i,j) == 1) then
1947 sum_vec_3d(i,j,1) = sum_vec_3d(i,j,1) + resid_scale * (zv(i,j) * rv(i,j))
1948 sum_vec_3d(i,j,2) = sum_vec_3d(i,j,2) + resid_scale * (dv(i,j) * av(i,j))
1949 rv_old(i,j) = rv(i,j) ; zv_old(i,j) = zv(i,j)
1950 endif
1951 enddo ; enddo
1952
1953 sv3dsum = reproducing_sum( sum_vec_3d(:,:,1:2), is_sum, ie_sum, js_sum, je_sum, sums=sv3dsums(1:2) )
1954
1955 alpha_k = sv3dsums(1)/sv3dsums(2)
1956
1957 do j=js,je-1 ; do i=is,ie-1
1958 if (cs%umask(i,j) == 1) then
1959 u_shlf(i,j) = u_shlf(i,j) + alpha_k * du(i,j)
1960 ru(i,j) = ru(i,j) - alpha_k * au(i,j)
1961 if (diagu(i,j)/=0) zu(i,j) = ru(i,j) / diagu(i,j)
1962 endif
1963 if (cs%vmask(i,j) == 1) then
1964 v_shlf(i,j) = v_shlf(i,j) + alpha_k * dv(i,j)
1965 rv(i,j) = rv(i,j) - alpha_k * av(i,j)
1966 if (diagv(i,j)/=0) zv(i,j) = rv(i,j) / diagv(i,j)
1967 endif
1968 enddo ; enddo
1969
1970
1971 ! R,u,v,Z valid region moves in by 1
1972
1973 ! beta_k = (Z \dot R) / (Zold \dot Rold)
1974 sum_vec_3d(:,:,:) = 0.0 ; sv3dsums(:) = 0.0
1975
1976 do j=jscq_sv,jecq ; do i=iscq_sv,iecq
1977 if (cs%umask(i,j) == 1) then
1978 sum_vec_3d(i,j,1) = resid_scale * (zu(i,j) * ru(i,j))
1979 sum_vec_3d(i,j,2) = resid_scale * (zu_old(i,j) * ru_old(i,j))
1980 sum_vec_3d(i,j,3) = resid2_scale * ru(i,j)**2
1981 endif
1982 if (cs%vmask(i,j) == 1) then
1983 sum_vec_3d(i,j,1) = sum_vec_3d(i,j,1) + resid_scale * (zv(i,j) * rv(i,j))
1984 sum_vec_3d(i,j,2) = sum_vec_3d(i,j,2) + resid_scale * (zv_old(i,j) * rv_old(i,j))
1985 sum_vec_3d(i,j,3) = sum_vec_3d(i,j,3) + resid2_scale * rv(i,j)**2
1986 endif
1987 enddo ; enddo
1988
1989 sv3dsum = reproducing_sum( sum_vec_3d, is_sum, ie_sum, js_sum, je_sum, sums=sv3dsums )
1990
1991 beta_k = sv3dsums(1)/sv3dsums(2)
1992
1993 do j=js,je-1 ; do i=is,ie-1
1994 if (cs%umask(i,j) == 1) du(i,j) = zu(i,j) + beta_k * du(i,j)
1995 if (cs%vmask(i,j) == 1) dv(i,j) = zv(i,j) + beta_k * dv(i,j)
1996 enddo ; enddo
1997
1998 ! D valid region moves in by 1
1999
2000 !if sqrt(sv3dsums(3)) <= (CS%cg_tolerance * resid0)
2001 if (sv3dsums(3) <= resid0tol2) then
2002 iters = iter
2003 conv_flag = 1
2004 exit
2005 endif
2006
2007 cg_halo = cg_halo - 1
2008
2009 if (cg_halo == 0) then
2010 ! pass vectors
2011 call pass_vector(du, dv, g%domain, to_all, bgrid_ne, complete=.false.)
2012 call pass_vector(u_shlf, v_shlf, g%domain, to_all, bgrid_ne, complete=.false.)
2013 call pass_vector(ru, rv, g%domain, to_all, bgrid_ne, complete=.true.)
2014 cg_halo = max_cg_halo
2015 endif
2016
2017 enddo ! end of CG loop
2018
2019 do j=jsdq,jedq ; do i=isdq,iedq
2020 if (cs%umask(i,j) == 3) then
2021 u_shlf(i,j) = cs%u_bdry_val(i,j)
2022 elseif (cs%umask(i,j) == 0) then
2023 u_shlf(i,j) = 0
2024 endif
2025
2026 if (cs%vmask(i,j) == 3) then
2027 v_shlf(i,j) = cs%v_bdry_val(i,j)
2028 elseif (cs%vmask(i,j) == 0) then
2029 v_shlf(i,j) = 0
2030 endif
2031 enddo ; enddo
2032
2033 call pass_vector(u_shlf, v_shlf, g%domain, to_all, bgrid_ne)
2034 if (conv_flag == 0) then
2035 iters = cs%cg_max_iterations
2036 endif
2037
2038end subroutine ice_shelf_solve_inner
2039
2040subroutine ice_shelf_advect_thickness_x(CS, G, LB, time_step, hmask, h0, h_after_uflux, uh_ice)
2041 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
2042 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
2043 type(loop_bounds_type), intent(in) :: LB !< Loop bounds structure.
2044 real, intent(in) :: time_step !< The time step for this update [T ~> s].
2045 real, dimension(SZDI_(G),SZDJ_(G)), &
2046 intent(inout) :: hmask !< A mask indicating which tracer points are
2047 !! partly or fully covered by an ice-shelf
2048 real, dimension(SZDI_(G),SZDJ_(G)), &
2049 intent(in) :: h0 !< The initial ice shelf thicknesses [Z ~> m].
2050 real, dimension(SZDI_(G),SZDJ_(G)), &
2051 intent(inout) :: h_after_uflux !< The ice shelf thicknesses after
2052 !! the zonal mass fluxes [Z ~> m].
2053 real, dimension(SZDIB_(G),SZDJ_(G)), &
2054 intent(inout) :: uh_ice !< The accumulated zonal ice volume flux [Z L2 ~> m3]
2055
2056 ! use will be made of ISS%hmask here - its value at the boundary will be zero, just like uncovered cells
2057 ! if there is an input bdry condition, the thickness there will be set in initialization
2058
2059
2060 integer :: i, j
2061 integer :: ish, ieh, jsh, jeh
2062 real :: u_face ! Zonal velocity at a face [L T-1 ~> m s-1]
2063 real :: h_face ! Thickness at a face for transport [Z ~> m]
2064 real :: slope_lim ! The value of the slope limiter, in the range of 0 to 2 [nondim]
2065
2066! is = G%isc-2 ; ie = G%iec+2 ; js = G%jsc ; je = G%jec
2067! isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
2068
2069 ish = lb%ish ; ieh = lb%ieh ; jsh = lb%jsh ; jeh = lb%jeh
2070
2071 ! hmask coded values: 1) fully covered; 2) partly covered - no export; 3) Specified boundary condition
2072 ! relevant u_face_mask coded values: 1) Normal interior point; 4) Specified flux BC
2073
2074 do j=jsh,jeh ; do i=ish-1,ieh
2075 if (cs%u_face_mask(i,j) == 4.) then ! The flux itself is a specified boundary condition.
2076 uh_ice(i,j) = (time_step * g%dyCu(i,j)) * cs%u_flux_bdry_val(i,j)
2077 elseif ((hmask(i,j) == 1 .or. hmask(i,j) == 3) .or. (hmask(i+1,j) == 1 .or. hmask(i+1,j) == 3)) then
2078 u_face = 0.5 * (cs%u_shelf(i,j-1) + cs%u_shelf(i,j))
2079 h_face = 0.0 ! This will apply when the source cell is iceless or not fully ice covered.
2080
2081 if (u_face > 0) then
2082 if (hmask(i,j) == 3) then ! This is a open boundary inflow from the west
2083 h_face = cs%h_bdry_val(i,j)
2084 elseif (hmask(i,j) == 1) then ! There can be eastward flow through this face.
2085 if ((hmask(i-1,j) == 1 .or. hmask(i-1,j) == 3) .and. &
2086 (hmask(i+1,j) == 1 .or. hmask(i+1,j) == 3)) then
2087 slope_lim = slope_limiter(h0(i,j)-h0(i-1,j), h0(i+1,j)-h0(i,j))
2088 ! This is a 2nd-order centered scheme with a slope limiter. We could try PPM here.
2089 h_face = h0(i,j) - slope_lim * (0.5 * (h0(i,j)-h0(i+1,j)))
2090 else
2091 h_face = h0(i,j)
2092 endif
2093 endif
2094 else
2095 if (hmask(i+1,j) == 3) then ! This is a open boundary inflow from the east
2096 h_face = cs%h_bdry_val(i+1,j)
2097 elseif (hmask(i+1,j) == 1) then
2098 if ((hmask(i,j) == 1 .or. hmask(i,j) == 3) .and. &
2099 (hmask(i+2,j) == 1 .or. hmask(i+2,j) == 3)) then
2100 slope_lim = slope_limiter(h0(i+1,j)-h0(i,j), h0(i+2,j)-h0(i+1,j))
2101 h_face = h0(i+1,j) - slope_lim * (0.5 * (h0(i+2,j)-h0(i+1,j)))
2102 else
2103 h_face = h0(i+1,j)
2104 endif
2105 endif
2106 endif
2107
2108 uh_ice(i,j) = (time_step * g%dyCu(i,j)) * (u_face * h_face)
2109 else
2110 uh_ice(i,j) = 0.0
2111 endif
2112 enddo ; enddo
2113
2114 do j=jsh,jeh ; do i=ish,ieh
2115 if (hmask(i,j) /= 3) &
2116 h_after_uflux(i,j) = h0(i,j) + (uh_ice(i-1,j) - uh_ice(i,j)) * g%IareaT(i,j)
2117
2118 ! Update the masks of cells that have gone from no ice to partial ice.
2119 if ((hmask(i,j) == 0) .and. ((uh_ice(i-1,j) > 0.0) .or. (uh_ice(i,j) < 0.0))) hmask(i,j) = 2
2120 enddo ; enddo
2121
2122end subroutine ice_shelf_advect_thickness_x
2123
2124subroutine ice_shelf_advect_thickness_y(CS, G, LB, time_step, hmask, h0, h_after_vflux, vh_ice)
2125 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
2126 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
2127 type(loop_bounds_type), intent(in) :: LB !< Loop bounds structure.
2128 real, intent(in) :: time_step !< The time step for this update [T ~> s].
2129 real, dimension(SZDI_(G),SZDJ_(G)), &
2130 intent(inout) :: hmask !< A mask indicating which tracer points are
2131 !! partly or fully covered by an ice-shelf
2132 real, dimension(SZDI_(G),SZDJ_(G)), &
2133 intent(in) :: h0 !< The initial ice shelf thicknesses [Z ~> m].
2134 real, dimension(SZDI_(G),SZDJ_(G)), &
2135 intent(inout) :: h_after_vflux !< The ice shelf thicknesses after
2136 !! the meridional mass fluxes [Z ~> m].
2137 real, dimension(SZDI_(G),SZDJB_(G)), &
2138 intent(inout) :: vh_ice !< The accumulated meridional ice volume flux [Z L2 ~> m3]
2139
2140 ! use will be made of ISS%hmask here - its value at the boundary will be zero, just like uncovered cells
2141 ! if there is an input bdry condition, the thickness there will be set in initialization
2142
2143
2144 integer :: i, j
2145 integer :: ish, ieh, jsh, jeh
2146 real :: v_face ! Pseudo-meridional velocity at a face [L T-1 ~> m s-1]
2147 real :: h_face ! Thickness at a face for transport [Z ~> m]
2148 real :: slope_lim ! The value of the slope limiter, in the range of 0 to 2 [nondim]
2149
2150 ish = lb%ish ; ieh = lb%ieh ; jsh = lb%jsh ; jeh = lb%jeh
2151
2152 ! hmask coded values: 1) fully covered; 2) partly covered - no export; 3) Specified boundary condition
2153 ! relevant u_face_mask coded values: 1) Normal interior point; 4) Specified flux BC
2154
2155 do j=jsh-1,jeh ; do i=ish,ieh
2156 if (cs%v_face_mask(i,j) == 4.) then ! The flux itself is a specified boundary condition.
2157 vh_ice(i,j) = (time_step * g%dxCv(i,j)) * cs%v_flux_bdry_val(i,j)
2158 elseif ((hmask(i,j) == 1 .or. hmask(i,j) == 3) .or. (hmask(i,j+1) == 1 .or. hmask(i,j+1) == 3)) then
2159 v_face = 0.5 * (cs%v_shelf(i-1,j) + cs%v_shelf(i,j))
2160 h_face = 0.0 ! This will apply when the source cell is iceless or not fully ice covered.
2161
2162 if (v_face > 0) then
2163 if (hmask(i,j) == 3) then ! This is a open boundary inflow from the south
2164 h_face = cs%h_bdry_val(i,j)
2165 elseif (hmask(i,j) == 1) then ! There can be northward flow through this face.
2166 if ((hmask(i,j-1) == 1 .or. hmask(i,j-1) == 3) .and. &
2167 (hmask(i,j+1) == 1 .or. hmask(i,j+1) == 3)) then
2168 slope_lim = slope_limiter(h0(i,j)-h0(i,j-1), h0(i,j+1)-h0(i,j))
2169 ! This is a 2nd-order centered scheme with a slope limiter. We could try PPM here.
2170 h_face = h0(i,j) - slope_lim * (0.5 * (h0(i,j)-h0(i,j+1)))
2171 else
2172 h_face = h0(i,j)
2173 endif
2174 endif
2175 else
2176 if (hmask(i,j+1) == 3) then ! This is a open boundary inflow from the north
2177 h_face = cs%h_bdry_val(i,j+1)
2178 elseif (hmask(i,j+1) == 1) then
2179 if ((hmask(i,j) == 1 .or. hmask(i,j) == 3) .and. &
2180 (hmask(i,j+2) == 1 .or. hmask(i,j+2) == 3)) then
2181 slope_lim = slope_limiter(h0(i,j+1)-h0(i,j), h0(i,j+2)-h0(i,j+1))
2182 h_face = h0(i,j+1) - slope_lim * (0.5 * (h0(i,j+2)-h0(i,j+1)))
2183 else
2184 h_face = h0(i,j+1)
2185 endif
2186 endif
2187 endif
2188
2189 vh_ice(i,j) = (time_step * g%dxCv(i,j)) * (v_face * h_face)
2190 else
2191 vh_ice(i,j) = 0.0
2192 endif
2193 enddo ; enddo
2194
2195 do j=jsh,jeh ; do i=ish,ieh
2196 if (hmask(i,j) /= 3) &
2197 h_after_vflux(i,j) = h0(i,j) + (vh_ice(i,j-1) - vh_ice(i,j)) * g%IareaT(i,j)
2198
2199 ! Update the masks of cells that have gone from no ice to partial ice.
2200 if ((hmask(i,j) == 0) .and. ((vh_ice(i,j-1) > 0.0) .or. (vh_ice(i,j) < 0.0))) hmask(i,j) = 2
2201 enddo ; enddo
2202
2203end subroutine ice_shelf_advect_thickness_y
2204
2205subroutine shelf_advance_front(CS, ISS, G, hmask, uh_ice, vh_ice)
2206 type(ice_shelf_dyn_cs), intent(in) :: cs !< A pointer to the ice shelf control structure
2207 type(ice_shelf_state), intent(inout) :: iss !< A structure with elements that describe
2208 !! the ice-shelf state
2209 type(ocean_grid_type), intent(in) :: g !< The grid structure used by the ice shelf.
2210 real, dimension(SZDI_(G),SZDJ_(G)), &
2211 intent(inout) :: hmask !< A mask indicating which tracer points are
2212 !! partly or fully covered by an ice-shelf
2213 real, dimension(SZDIB_(G),SZDJ_(G)), &
2214 intent(inout) :: uh_ice !< The accumulated zonal ice volume flux [Z L2 ~> m3]
2215 real, dimension(SZDI_(G),SZDJB_(G)), &
2216 intent(inout) :: vh_ice !< The accumulated meridional ice volume flux [Z L2 ~> m3]
2217
2218 ! in this subroutine we go through the computational cells only and, if they are empty or partial cells,
2219 ! we find the reference thickness and update the shelf mass and partial area fraction and the hmask if necessary
2220
2221 ! if any cells go from partial to complete, we then must set the thickness, update hmask accordingly,
2222 ! and divide the overflow across the adjacent EMPTY (not partly-covered) cells.
2223 ! (it is highly unlikely there will not be any; in which case this will need to be rethought.)
2224
2225 ! most likely there will only be one "overflow". If not, though, a pass_var of all relevant variables
2226 ! is done; there will therefore be a loop which, in practice, will hopefully not have to go through
2227 ! many iterations
2228
2229 ! when 3d advected scalars are introduced, they will be impacted by what is done here
2230
2231 ! flux_enter(isd:ied,jsd:jed,1:4): if cell is not ice-covered, gives flux of ice into cell from kth boundary
2232 !
2233 ! from eastern neighbor: flux_enter(:,:,1)
2234 ! from western neighbor: flux_enter(:,:,2)
2235 ! from southern neighbor: flux_enter(:,:,3)
2236 ! from northern neighbor: flux_enter(:,:,4)
2237 !
2238 ! o--- (4) ---o
2239 ! | |
2240 ! (1) (2)
2241 ! | |
2242 ! o--- (3) ---o
2243 !
2244
2245 integer :: i, j, isc, iec, jsc, jec, n_flux, k, iter_count
2246 integer :: i_off, j_off
2247 integer :: iter_flag
2248
2249 real :: h_reference ! A reference thicknesss based on neighboring cells [Z ~> m]
2250 real :: h_reference_ew !contribution to reference thickness from east + west cells [Z ~> m]
2251 real :: h_reference_ns !contribution to reference thickness from north + south cells [Z ~> m]
2252 real :: tot_flux ! The total ice mass flux [Z L2 ~> m3]
2253 real :: tot_flux_ew ! The contribution to total ice mass flux from east + west cells [Z L2 ~> m3]
2254 real :: tot_flux_ns ! The contribution to total ice mass flux from north + south cells [Z L2 ~> m3]
2255 real :: partial_vol ! The volume covered by ice shelf [Z L2 ~> m3]
2256 real :: dxdyh ! Cell area [L2 ~> m2]
2257 character(len=160) :: mesg ! The text of an error message
2258 integer, dimension(4) :: mapi, mapj, new_partial
2259 real, dimension(SZDI_(G),SZDJ_(G),4) :: flux_enter ! The ice volume flux into the
2260 ! cell through the 4 cell boundaries [Z L2 ~> m3].
2261 real, dimension(SZDI_(G),SZDJ_(G),4) :: flux_enter_replace ! An updated ice volume flux into the
2262 ! cell through the 4 cell boundaries [Z L2 ~> m3].
2263
2264 isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
2265 i_off = g%idg_offset ; j_off = g%jdg_offset
2266 iter_count = 0 ; iter_flag = 1
2267
2268 flux_enter(:,:,:) = 0.0
2269 do j=jsc-1,jec+1 ; do i=isc-1,iec+1
2270 if ((hmask(i,j) == 0) .or. (hmask(i,j) == 2)) then
2271 flux_enter(i,j,1) = max(uh_ice(i-1,j), 0.0)
2272 flux_enter(i,j,2) = max(-uh_ice(i,j), 0.0)
2273 flux_enter(i,j,3) = max(vh_ice(i,j-1), 0.0)
2274 flux_enter(i,j,4) = max(-vh_ice(i,j), 0.0)
2275 endif
2276 enddo ; enddo
2277
2278 mapi(1) = -1 ; mapi(2) = 1 ; mapi(3:4) = 0
2279 mapj(3) = -1 ; mapj(4) = 1 ; mapj(1:2) = 0
2280
2281 do while (iter_flag == 1)
2282
2283 iter_flag = 0
2284
2285 if (iter_count > 0) then
2286 flux_enter(:,:,:) = flux_enter_replace(:,:,:)
2287 endif
2288 flux_enter_replace(:,:,:) = 0.0
2289
2290 iter_count = iter_count + 1
2291
2292 ! if iter_count >= 3 then some halo updates need to be done...
2293 if (iter_count==3) then
2294 call mom_error(fatal, "MOM_ice_shelf_dyn.F90, shelf_advance_front iter >=3.")
2295 endif
2296
2297 do j=jsc-1,jec+1
2298
2299 if (cs%reentrant_y .OR. (((j+j_off) <= g%domain%njglobal) .AND. &
2300 ((j+j_off) >= 1))) then
2301
2302 do i=isc-1,iec+1
2303
2304 if (cs%reentrant_x .OR. (((i+i_off) <= g%domain%niglobal) .AND. &
2305 ((i+i_off) >= 1))) then
2306 ! first get reference thickness by averaging over cells that are fluxing into this cell
2307 n_flux = 0
2308 h_reference_ew = 0.0
2309 h_reference_ns = 0.0
2310 tot_flux_ew = 0.0
2311 tot_flux_ns = 0.0
2312
2313 do k=1,2
2314 if (flux_enter(i,j,k) > 0) then
2315 n_flux = n_flux + 1
2316 h_reference_ew = h_reference_ew + flux_enter(i,j,k) * iss%h_shelf(i+2*k-3,j)
2317 !h_reference = h_reference + ISS%h_shelf(i+2*k-3,j)
2318 tot_flux_ew = tot_flux_ew + flux_enter(i,j,k)
2319 flux_enter(i,j,k) = 0.0
2320 endif
2321 enddo
2322
2323 do k=1,2
2324 if (flux_enter(i,j,k+2) > 0) then
2325 n_flux = n_flux + 1
2326 h_reference_ns = h_reference_ns + flux_enter(i,j,k+2) * iss%h_shelf(i,j+2*k-3)
2327 !h_reference = h_reference + ISS%h_shelf(i,j+2*k-3)
2328 tot_flux_ns = tot_flux_ns + flux_enter(i,j,k+2)
2329 flux_enter(i,j,k+2) = 0.0
2330 endif
2331 enddo
2332
2333 h_reference = h_reference_ew + h_reference_ns
2334 tot_flux = tot_flux_ew + tot_flux_ns
2335
2336 if (n_flux > 0) then
2337 dxdyh = g%areaT(i,j)
2338 h_reference = h_reference / tot_flux
2339 !h_reference = h_reference / real(n_flux)
2340 partial_vol = iss%h_shelf(i,j) * iss%area_shelf_h(i,j) + tot_flux
2341
2342 if ((partial_vol / g%areaT(i,j)) == h_reference) then ! cell is exactly covered, no overflow
2343 if (iss%hmask(i,j)/=3) iss%hmask(i,j) = 1
2344 iss%h_shelf(i,j) = h_reference
2345 iss%area_shelf_h(i,j) = g%areaT(i,j)
2346 elseif ((partial_vol / g%areaT(i,j)) < h_reference) then
2347 iss%hmask(i,j) = 2
2348 ! ISS%mass_shelf(i,j) = partial_vol * CS%density_ice
2349 iss%area_shelf_h(i,j) = partial_vol / h_reference
2350 iss%h_shelf(i,j) = h_reference
2351 else
2352
2353 if (iss%hmask(i,j)/=3) iss%hmask(i,j) = 1
2354 iss%area_shelf_h(i,j) = g%areaT(i,j)
2355 !h_temp(i,j) = h_reference
2356 partial_vol = partial_vol - h_reference * g%areaT(i,j)
2357
2358 iter_flag = 1
2359
2360 n_flux = 0 ; new_partial(:) = 0
2361
2362 do k=1,2
2363 if (cs%u_face_mask(i-2+k,j) == 2) then
2364 n_flux = n_flux + 1
2365 elseif (iss%hmask(i+2*k-3,j) == 0) then
2366 n_flux = n_flux + 1
2367 new_partial(k) = 1
2368 endif
2369 if (cs%v_face_mask(i,j-2+k) == 2) then
2370 n_flux = n_flux + 1
2371 elseif (iss%hmask(i,j+2*k-3) == 0) then
2372 n_flux = n_flux + 1
2373 new_partial(k+2) = 1
2374 endif
2375 enddo
2376
2377 if (n_flux == 0) then ! there is nowhere to put the extra ice!
2378 iss%h_shelf(i,j) = h_reference + partial_vol / g%areaT(i,j)
2379 else
2380 iss%h_shelf(i,j) = h_reference
2381
2382 do k=1,2
2383 if (new_partial(k) == 1) &
2384 flux_enter_replace(i+2*k-3,j,3-k) = partial_vol / real(n_flux)
2385 if (new_partial(k+2) == 1) &
2386 flux_enter_replace(i,j+2*k-3,5-k) = partial_vol / real(n_flux)
2387 enddo
2388 endif
2389
2390 endif ! Parital_vol test.
2391 endif ! n_flux gt 0 test.
2392
2393 endif
2394 enddo ! j-loop
2395 endif
2396 enddo
2397
2398 ! call max_across_PEs(iter_flag)
2399
2400 enddo ! End of do while(iter_flag) loop
2401
2402 call max_across_pes(iter_count)
2403
2404 if (is_root_pe() .and. (iter_count > 1)) then
2405 write(mesg,*) "shelf_advance_front: ", iter_count, " max iterations"
2406 call mom_mesg(mesg, 5)
2407 endif
2408
2409end subroutine shelf_advance_front
2410
2411!> Apply a very simple calving law using a minimum thickness rule
2412subroutine ice_shelf_min_thickness_calve(G, h_shelf, area_shelf_h, hmask, thickness_calve, halo)
2413 type(ocean_grid_type), intent(in) :: g !< The grid structure used by the ice shelf.
2414 real, dimension(SZDI_(G),SZDJ_(G)), intent(inout) :: h_shelf !< The ice shelf thickness [Z ~> m].
2415 real, dimension(SZDI_(G),SZDJ_(G)), intent(inout) :: area_shelf_h !< The area per cell covered by
2416 !! the ice shelf [L2 ~> m2].
2417 real, dimension(SZDI_(G),SZDJ_(G)), intent(inout) :: hmask !< A mask indicating which tracer points are
2418 !! partly or fully covered by an ice-shelf
2419 real, intent(in) :: thickness_calve !< The thickness at which to trigger calving [Z ~> m].
2420 integer, optional, intent(in) :: halo !< The number of halo points to use. If not present,
2421 !! work on the entire data domain.
2422 integer :: i, j, is, ie, js, je
2423
2424 if (present(halo)) then
2425 is = g%isc - halo ; ie = g%iec + halo ; js = g%jsc - halo ; je = g%jec + halo
2426 else
2427 is = g%isd ; ie = g%ied ; js = g%jsd ; je = g%jed
2428 endif
2429
2430 do j=js,je ; do i=is,ie
2431! if ((h_shelf(i,j) < CS%thickness_calve) .and. (hmask(i,j) == 1) .and. &
2432! (CS%ground_frac(i,j) == 0.0)) then
2433 if ((h_shelf(i,j) < thickness_calve) .and. (area_shelf_h(i,j) > 0.)) then
2434 h_shelf(i,j) = 0.0
2435 area_shelf_h(i,j) = 0.0
2436 hmask(i,j) = 0.0
2437 endif
2438 enddo ; enddo
2439
2440end subroutine ice_shelf_min_thickness_calve
2441
2442subroutine calve_to_mask(G, h_shelf, area_shelf_h, hmask, calve_mask)
2443 type(ocean_grid_type), intent(in) :: g !< The grid structure used by the ice shelf.
2444 real, dimension(SZDI_(G),SZDJ_(G)), intent(inout) :: h_shelf !< The ice shelf thickness [Z ~> m].
2445 real, dimension(SZDI_(G),SZDJ_(G)), intent(inout) :: area_shelf_h !< The area per cell covered by
2446 !! the ice shelf [L2 ~> m2].
2447 real, dimension(SZDI_(G),SZDJ_(G)), intent(inout) :: hmask !< A mask indicating which tracer points are
2448 !! partly or fully covered by an ice-shelf
2449 real, dimension(SZDI_(G),SZDJ_(G)), intent(in) :: calve_mask !< A mask that indicates where the ice
2450 !! shelf can exist, and where it will calve.
2451
2452 integer :: i,j
2453
2454 do j=g%jsc,g%jec ; do i=g%isc,g%iec
2455 if ((calve_mask(i,j) == 0.0) .and. (hmask(i,j) /= 0.0)) then
2456 h_shelf(i,j) = 0.0
2457 area_shelf_h(i,j) = 0.0
2458 hmask(i,j) = 0.0
2459 endif
2460 enddo ; enddo
2461
2462end subroutine calve_to_mask
2463
2464!> Calculate driving stress using cell-centered bed elevation and ice thickness
2465subroutine calc_shelf_driving_stress(CS, ISS, G, US, taudx, taudy, OD)
2466 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
2467 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
2468 !! the ice-shelf state
2469 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
2470 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
2471 real, dimension(SZDI_(G),SZDJ_(G)), &
2472 intent(in) :: OD !< ocean floor depth at tracer points [Z ~> m].
2473 real, dimension(SZDIB_(G),SZDJB_(G)), &
2474 intent(inout) :: taudx !< X-direction driving stress at q-points [R L3 Z T-2 ~> kg m s-2]
2475 real, dimension(SZDIB_(G),SZDJB_(G)), &
2476 intent(inout) :: taudy !< Y-direction driving stress at q-points [R L3 Z T-2 ~> kg m s-2]
2477
2478
2479! driving stress!
2480
2481! ! taudx and taudy will hold driving stress in the x- and y- directions when done.
2482! they will sit on the BGrid, and so their size depends on whether the grid is symmetric
2483!
2484! Since this is a finite element solve, they will actually have the form \int \Phi_i rho g h \nabla s
2485!
2486! OD -this is important and we do not yet know where (in MOM) it will come from. It represents
2487! "average" ocean depth -- and is needed to find surface elevation
2488! (it is assumed that base_ice = bed + OD)
2489
2490 real, dimension(SIZE(OD,1),SIZE(OD,2)) :: S ! surface elevation [Z ~> m].
2491 real, dimension(SZDI_(G),SZDJ_(G)) :: sx_e, sy_e !element contributions to driving stress
2492 real :: rho, rhow, rhoi_rhow ! Ice and ocean densities [R ~> kg m-3]
2493 real :: sx, sy ! Ice shelf top slopes at tracer points [Z L-1 ~> nondim]
2494 real :: neumann_val ! [R Z L2 T-2 ~> kg s-2]
2495 real :: grav ! The gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
2496 real :: scale ! Scaling factor used to ensure surface slope magnitude does not exceed CS%max_surface_slope
2497 logical :: valid_N, valid_S, valid_E, valid_W
2498 integer :: i, j, iscq, iecq, jscq, jecq, isd, jsd, ied, jed, is, js, iegq, jegq
2499 integer :: giec, gjec, gisc, gjsc, isc, jsc, iec, jec
2500 integer :: i_off, j_off
2501
2502 isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
2503! iscq = G%iscB ; iecq = G%iecB ; jscq = G%jscB ; jecq = G%jecB
2504 isd = g%isd ; jsd = g%jsd ; ied = g%ied ; jed = g%jed
2505! iegq = G%iegB ; jegq = G%jegB
2506! gisc = G%domain%nihalo+1 ; gjsc = G%domain%njhalo+1
2507 gisc = 1 ; gjsc = 1
2508! giec = G%domain%niglobal+G%domain%nihalo ; gjec = G%domain%njglobal+G%domain%njhalo
2509 giec = g%domain%niglobal ; gjec = g%domain%njglobal
2510! is = iscq - 1 ; js = jscq - 1
2511 i_off = g%idg_offset ; j_off = g%jdg_offset
2512
2513
2514 rho = cs%density_ice
2515 rhow = cs%density_ocean_avg
2516 grav = cs%g_Earth
2517 rhoi_rhow = rho/rhow
2518 ! prelim - go through and calculate S
2519
2520 if (cs%GL_couple) then
2521 do j=jsc-2,jec+2 ; do i=isc-2,iec+2
2522 s(i,j) = -cs%bed_elev(i,j) + (od(i,j) + max(iss%h_shelf(i,j),cs%min_h_shelf))
2523 enddo ; enddo
2524 else
2525 ! check whether the ice is floating or grounded
2526 do j=jsc-2,jec+2 ; do i=isc-2,iec+2
2527 if (rhoi_rhow * max(iss%h_shelf(i,j),cs%min_h_shelf) - cs%bed_elev(i,j) <= 0) then
2528 s(i,j) = (1 - rhoi_rhow)*max(iss%h_shelf(i,j),cs%min_h_shelf)
2529 else
2530 s(i,j) = max(iss%h_shelf(i,j),cs%min_h_shelf)-cs%bed_elev(i,j)
2531 endif
2532 enddo ; enddo
2533 endif
2534
2535 call pass_var(s, g%domain)
2536
2537 do j=jsc-1,jec+1
2538 do i=isc-1,iec+1
2539
2540 if (iss%hmask(i,j) == 1 .or. iss%hmask(i,j) == 3) then
2541 ! we are inside the global computational bdry, at an ice-filled cell
2542
2543 ! Calculate the x-direction surface slope at tracer points.
2544 sx = 0.0
2545 valid_e = (iss%hmask(i+1,j) == 1 .or. iss%hmask(i+1,j) == 3)
2546 valid_w = (iss%hmask(i-1,j) == 1 .or. iss%hmask(i-1,j) == 3)
2547 if (cs%shelf_top_slope_bugs) then
2548 if (((i+i_off) == gisc) .and. (.not.cs%reentrant_x)) then ! at west computational bdry
2549 if (valid_e) sx = (s(i+1,j)-s(i,j)) / g%dxT(i,j)
2550 elseif (((i+i_off) == giec) .and. (.not.cs%reentrant_x)) then ! at east computational bdry
2551 if (valid_w) sx = (s(i,j)-s(i-1,j)) / g%dxT(i,j)
2552 elseif (valid_e .and. valid_w) then
2553 ! This is the usual interior point
2554 sx = (s(i+1,j) - s(i-1,j)) / (g%dxT(i,j) + g%dxT(i-1,j))
2555 elseif (valid_e) then
2556 sx = (s(i+1,j) - s(i,j)) / (g%dxT(i,j) + g%dxT(i+1,j))
2557 elseif (valid_w) then
2558 sx = (s(i,j) - s(i-1,j)) / (g%dxT(i,j) + g%dxT(i-1,j))
2559 endif
2560 else ! Correct the bugs in the version above.
2561 if (((i+i_off) == gisc) .and. (.not.cs%reentrant_x)) then ! at west computational bdry
2562 if (valid_e) sx = (s(i+1,j) - s(i,j)) * g%IdxCu(i,j)
2563 elseif (((i+i_off) == giec) .and. (.not.cs%reentrant_x)) then ! at east computational bdry
2564 if (valid_w) sx = (s(i,j) - s(i-1,j)) * g%IdxCu(i-1,j)
2565 elseif (valid_e .and. valid_w) then
2566 ! This is the usual interior point
2567 sx = 0.5*(s(i+1,j) - s(i-1,j)) * g%IdxT(i,j)
2568 elseif (valid_e) then ! Use a one-sided estimate from the east.
2569 sx = (s(i+1,j) - s(i,j)) * g%IdxCu(i,j)
2570 elseif (valid_w) then ! Use a one-sided estimate from the west.
2571 sx = (s(i,j) - s(i-1,j)) * g%IdxCu(i-1,j)
2572 endif
2573 endif
2574
2575 ! Calculate the y-direction surface slope at tracer points.
2576 sy = 0.0
2577 valid_n = (iss%hmask(i,j+1) == 1 .or. iss%hmask(i,j+1) == 3)
2578 valid_s = (iss%hmask(i,j-1) == 1 .or. iss%hmask(i,j-1) == 3)
2579 if (cs%shelf_top_slope_bugs) then
2580 if (((j+j_off) == gjsc) .and. (.not. cs%reentrant_y)) then ! at south computational bdry
2581 if (valid_n) sy = (s(i,j+1)-s(i,j)) / g%dyT(i,j)
2582 elseif (((j+j_off) == gjec) .and. (.not. cs%reentrant_y)) then ! at north computational bdry
2583 if (valid_s) sy = (s(i,j)-s(i,j-1)) / g%dyT(i,j)
2584 elseif (valid_n .and. valid_s) then
2585 ! This is the usual interior point
2586 sy = (s(i,j+1) - s(i,j-1)) / (g%dyT(i,j) + g%dyT(i,j-1))
2587 elseif (valid_n) then
2588 sy = (s(i,j+1) - s(i,j)) / (g%dyT(i,j) + g%dyT(i,j+1))
2589 elseif (valid_s) then
2590 sy = (s(i,j) - s(i,j-1)) / (g%dyT(i,j) + g%dyT(i,j-1))
2591 endif
2592 else ! Correct the bugs in the version above.
2593 if (((j+j_off) == gjsc) .and. (.not. cs%reentrant_y)) then ! at south computational bdry
2594 if (valid_n) sy = (s(i,j+1) - s(i,j)) * g%IdyCv(i,j)
2595 elseif (((j+j_off) == gjec) .and. (.not. cs%reentrant_y)) then ! at north computational bdry
2596 if (valid_s) sy = (s(i,j) - s(i,j-1)) * g%IdyCv(i,j-1)
2597 elseif (valid_n .and. valid_s) then
2598 ! This is the usual interior point
2599 sy = 0.5*(s(i,j+1) - s(i,j-1)) * g%IdyT(i,j)
2600 elseif (valid_n) then ! Use a one-sided estimate from the north.
2601 sy = (s(i,j+1) - s(i,j)) * g%IdyCv(i,j)
2602 elseif (valid_s) then ! Use a one-sided estimate from the south.
2603 sy = (s(i,j) - s(i,j-1)) * g%IdyCv(i,j-1)
2604 endif
2605 endif
2606
2607 if (cs%max_surface_slope>0) then
2608 scale = cs%max_surface_slope / max( sqrt((sx**2) + (sy**2)), cs%max_surface_slope )
2609 sx = scale*sx ; sy = scale*sy
2610 endif
2611
2612 sx_e(i,j) = (-.25 * g%areaT(i,j)) * ((rho * grav) * (max(iss%h_shelf(i,j),cs%min_h_shelf) * sx))
2613 sy_e(i,j) = (-.25 * g%areaT(i,j)) * ((rho * grav) * (max(iss%h_shelf(i,j),cs%min_h_shelf) * sy))
2614
2615 cs%sx_shelf(i,j) = sx ; cs%sy_shelf(i,j) = sy
2616
2617 !Stress (Neumann) boundary conditions
2618 if (cs%ground_frac(i,j) == 1) then
2619 neumann_val = ((.5 * grav) * (rho * max(iss%h_shelf(i,j),cs%min_h_shelf)**2 - rhow * cs%bed_elev(i,j)**2))
2620 else
2621 neumann_val = (.5 * grav) * ((1-rho/rhow) * (rho * max(iss%h_shelf(i,j),cs%min_h_shelf)**2))
2622 endif
2623 if ((cs%u_face_mask_bdry(i-1,j) == 2) .OR. &
2624 ((iss%hmask(i-1,j) == 0 .OR. iss%hmask(i-1,j) == 2) .AND. (cs%reentrant_x .OR. (i+i_off /= gisc)))) then
2625 ! left face of the cell is at a stress boundary
2626 ! the depth-integrated longitudinal stress is equal to the difference of depth-integrated
2627 ! pressure on either side of the face
2628 ! on the ice side, it is rho g h^2 / 2
2629 ! on the ocean side, it is rhow g (delta OD)^2 / 2
2630 ! OD can be zero under the ice; but it is ASSUMED on the ice-free side of the face, topography elevation
2631 ! is not above the base of the ice in the current cell
2632
2633 ! Note the negative sign due to the direction of the normal vector
2634 taudx(i-1,j-1) = taudx(i-1,j-1) - .5 * g%dyT(i,j) * neumann_val
2635 taudx(i-1,j) = taudx(i-1,j) - .5 * g%dyT(i,j) * neumann_val
2636 endif
2637
2638 if ((cs%u_face_mask_bdry(i,j) == 2) .OR. &
2639 ((iss%hmask(i+1,j) == 0 .OR. iss%hmask(i+1,j) == 2) .and. (cs%reentrant_x .OR. (i+i_off /= giec)))) then
2640 ! east face of the cell is at a stress boundary
2641 taudx(i,j-1) = taudx(i,j-1) + .5 * g%dyT(i,j) * neumann_val
2642 taudx(i,j) = taudx(i,j) + .5 * g%dyT(i,j) * neumann_val
2643 endif
2644
2645 if ((cs%v_face_mask_bdry(i,j-1) == 2) .OR. &
2646 ((iss%hmask(i,j-1) == 0 .OR. iss%hmask(i,j-1) == 2) .and. (cs%reentrant_y .OR. (j+j_off /= gjsc)))) then
2647 ! south face of the cell is at a stress boundary
2648 taudy(i-1,j-1) = taudy(i-1,j-1) - .5 * g%dxT(i,j) * neumann_val
2649 taudy(i,j-1) = taudy(i,j-1) - .5 * g%dxT(i,j) * neumann_val
2650 endif
2651
2652 if ((cs%v_face_mask_bdry(i,j) == 2) .OR. &
2653 ((iss%hmask(i,j+1) == 0 .OR. iss%hmask(i,j+1) == 2) .and. (cs%reentrant_y .OR. (j+j_off /= gjec)))) then
2654 ! north face of the cell is at a stress boundary
2655 taudy(i-1,j) = taudy(i-1,j) + .5 * g%dxT(i,j) * neumann_val
2656 taudy(i,j) = taudy(i,j) + .5 * g%dxT(i,j) * neumann_val
2657 endif
2658 else ! This is not an ice-filled cell, so zero out the slopes here
2659 cs%sx_shelf(i,j) = 0.0 ; cs%sy_shelf(i,j) = 0.0
2660 sx_e(i,j) = 0.0
2661 sy_e(i,j) = 0.0
2662 endif
2663 enddo
2664 enddo
2665
2666 do j=jsc-1,jec ; do i=isc-1,iec
2667 taudx(i,j) = taudx(i,j) + ((sx_e(i,j)+sx_e(i+1,j+1)) + (sx_e(i+1,j)+sx_e(i,j+1)))
2668 taudy(i,j) = taudy(i,j) + ((sy_e(i,j)+sy_e(i+1,j+1)) + (sy_e(i+1,j)+sy_e(i,j+1)))
2669 enddo ; enddo
2670end subroutine calc_shelf_driving_stress
2671
2672subroutine cg_action(CS, uret, vret, u_shlf, v_shlf, Phi, Phisub, umask, vmask, hmask, H_node, &
2673 ice_visc, float_cond, bathyT, basal_trac, G, US, is, ie, js, je, dens_ratio, use_newton_in)
2674
2675 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
2676 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
2677 real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
2678 intent(inout) :: uret !< The retarding stresses working at u-points [R L3 Z T-2 ~> kg m s-2].
2679 real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
2680 intent(inout) :: vret !< The retarding stresses working at v-points [R L3 Z T-2 ~> kg m s-2].
2681 real, dimension(8,4,SZDI_(G),SZDJ_(G)), &
2682 intent(in) :: Phi !< The gradients of bilinear basis elements at Gaussian
2683 !! quadrature points surrounding the cell vertices [L-1 ~> m-1].
2684 real, dimension(:,:,:,:,:,:), &
2685 intent(in) :: Phisub !< Quadrature structure weights at subgridscale
2686 !! locations for finite element calculations [nondim]
2687 real, dimension(SZDIB_(G),SZDJB_(G)), &
2688 intent(in) :: u_shlf !< The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
2689 real, dimension(SZDIB_(G),SZDJB_(G)), &
2690 intent(in) :: v_shlf !< The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
2691 real, dimension(SZDIB_(G),SZDJB_(G)), &
2692 intent(in) :: umask !< A coded mask indicating the nature of the
2693 !! zonal flow at the corner point
2694 real, dimension(SZDIB_(G),SZDJB_(G)), &
2695 intent(in) :: vmask !< A coded mask indicating the nature of the
2696 !! meridional flow at the corner point
2697 real, dimension(SZDIB_(G),SZDJB_(G)), &
2698 intent(in) :: H_node !< The ice shelf thickness at nodal (corner)
2699 !! points [Z ~> m].
2700 real, dimension(SZDI_(G),SZDJ_(G)), &
2701 intent(in) :: hmask !< A mask indicating which tracer points are
2702 !! partly or fully covered by an ice-shelf
2703 real, dimension(SZDI_(G),SZDJ_(G),CS%visc_qps), &
2704 intent(in) :: ice_visc !< A field related to the ice viscosity from Glen's
2705 !! flow law [R L4 Z T-1 ~> kg m2 s-1].
2706 real, dimension(SZDI_(G),SZDJ_(G)), &
2707 intent(in) :: float_cond !< If GL_regularize=true, an array indicating where the ice
2708 !! shelf is floating: 0 if floating, 1 if not
2709 real, dimension(SZDI_(G),SZDJ_(G)), &
2710 intent(in) :: bathyT !< The depth of ocean bathymetry at tracer points
2711 !! relative to sea-level [Z ~> m].
2712 real, dimension(SZDI_(G),SZDJ_(G)), &
2713 intent(in) :: basal_trac !< Area-integrated taub_beta field related to the nonlinear
2714 !! part of the "linearized" basal stress [R Z L2 T-1 ~> kg s-1].
2715
2716 real, intent(in) :: dens_ratio !< The density of ice divided by the density
2717 !! of seawater, nondimensional
2718 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
2719 integer, intent(in) :: is !< The starting i-index to work on
2720 integer, intent(in) :: ie !< The ending i-index to work on
2721 integer, intent(in) :: js !< The starting j-index to work on
2722 integer, intent(in) :: je !< The ending j-index to work on
2723 logical, optional, intent(in) :: use_newton_in !< If present, overrides CS%doing_newton for Newton correction
2724
2725! the linear action of the matrix on (u,v) with bilinear finite elements
2726! as of now everything is passed in so no grid pointers or anything of the sort have to be dereferenced,
2727! but this may change pursuant to conversations with others
2728!
2729! is & ie are the cells over which the iteration is done; this may change between calls to this subroutine
2730! in order to make less frequent halo updates
2731
2732! the linear action of the matrix on (u,v) with bilinear finite elements
2733! Phi has the form
2734! Phi(k,q,i,j) - applies to cell i,j
2735
2736 ! 3 - 4
2737 ! | |
2738 ! 1 - 2
2739
2740! Phi(2*k-1,q,i,j) gives d(Phi_k)/dx at quadrature point q
2741! Phi(2*k,q,i,j) gives d(Phi_k)/dy at quadrature point q
2742! Phi_k is equal to 1 at vertex k, and 0 at vertex l /= k, and bilinear
2743
2744 real :: ux, uy, vx, vy ! Components of velocity shears or divergence [T-1 ~> s-1]
2745 real :: uq, vq ! Interpolated velocities [L T-1 ~> m s-1]
2746 real :: strx_n, stry_n, strsh_n, dstrain_n, inner_dot_n ! Newton correction variables [T-1 ~> s-1], [T-2 ~> s-2]
2747 integer :: iq, jq, iphi, jphi, i, j, ilq, jlq, Itgt, Jtgt, qp, qpv
2748 logical :: visc_qp4
2749 logical :: use_newton ! Whether to apply Newton tangent stiffness corrections
2750 logical :: do_newton_visc ! Whether to apply viscosity-related Newton tangent stiffness corrections
2751 real, dimension(2) :: xquad ! Nondimensional quadrature ratios [nondim]
2752 real, dimension(2,2) :: Ucell, Vcell, Usub, Vsub ! Velocities at the nodal points around the cell [L T-1 ~> m s-1]
2753 real, dimension(2,2) :: Hcell ! Ice shelf thickness at notal (corner) points [Z ~> m]
2754 real, dimension(2,2,4) :: uret_qp, vret_qp ! Temporary arrays in [R Z L3 T-2 ~> kg m s-2]
2755 real, dimension(SZDIB_(G),SZDJB_(G),4) :: uret_b, vret_b ! Temporary arrays in [R Z L3 T-2 ~> kg m s-2]
2756
2757 xquad(1) = .5 * (1-sqrt(1./3)) ; xquad(2) = .5 * (1+sqrt(1./3))
2758
2759 if (cs%visc_qps == 4) then
2760 visc_qp4=.true.
2761 else
2762 visc_qp4=.false.
2763 qpv = 1
2764 endif
2765
2766 use_newton = cs%doing_newton
2767 if (present(use_newton_in)) use_newton = use_newton_in
2768 do_newton_visc = use_newton .and. trim(cs%ice_viscosity_compute) == "MODEL"
2769
2770 uret(:,:) = 0.0 ; vret(:,:) = 0.0
2771 uret_b(:,:,:) = 0.0 ; vret_b(:,:,:) = 0.0
2772
2773 do j=js,je ; do i=is,ie ; if (hmask(i,j) == 1 .or. hmask(i,j)==3) then
2774
2775 uret_qp(:,:,:) = 0.0 ; vret_qp(:,:,:) = 0.0
2776
2777 do iq=1,2 ; do jq=1,2
2778
2779 qp = 2*(jq-1)+iq !current quad point
2780
2781 uq = ((u_shlf(i-1,j-1) * (xquad(3-iq) * xquad(3-jq))) + &
2782 (u_shlf(i,j) * (xquad(iq) * xquad(jq)))) + &
2783 ((u_shlf(i,j-1) * (xquad(iq) * xquad(3-jq))) + &
2784 (u_shlf(i-1,j) * (xquad(3-iq) * xquad(jq))))
2785
2786 vq = ((v_shlf(i-1,j-1) * (xquad(3-iq) * xquad(3-jq))) + &
2787 (v_shlf(i,j) * (xquad(iq) * xquad(jq)))) + &
2788 ((v_shlf(i,j-1) * (xquad(iq) * xquad(3-jq))) + &
2789 (v_shlf(i-1,j) * (xquad(3-iq) * xquad(jq))))
2790
2791 ux = ((u_shlf(i-1,j-1) * phi(1,qp,i,j)) + &
2792 (u_shlf(i,j) * phi(7,qp,i,j))) + &
2793 ((u_shlf(i,j-1) * phi(3,qp,i,j)) + &
2794 (u_shlf(i-1,j) * phi(5,qp,i,j)))
2795
2796 vx = ((v_shlf(i-1,j-1) * phi(1,qp,i,j)) + &
2797 (v_shlf(i,j) * phi(7,qp,i,j))) + &
2798 ((v_shlf(i,j-1) * phi(3,qp,i,j)) + &
2799 (v_shlf(i-1,j) * phi(5,qp,i,j)))
2800
2801 uy = ((u_shlf(i-1,j-1) * phi(2,qp,i,j)) + &
2802 (u_shlf(i,j) * phi(8,qp,i,j))) + &
2803 ((u_shlf(i,j-1) * phi(4,qp,i,j)) + &
2804 (u_shlf(i-1,j) * phi(6,qp,i,j)))
2805
2806 vy = ((v_shlf(i-1,j-1) * phi(2,qp,i,j)) + &
2807 (v_shlf(i,j) * phi(8,qp,i,j))) + &
2808 ((v_shlf(i,j-1) * phi(4,qp,i,j)) + &
2809 (v_shlf(i-1,j) * phi(6,qp,i,j)))
2810
2811 if (visc_qp4) qpv = qp !current quad point for viscosity
2812
2813 ! Newton correction: compute dstrain scalar once per quadrature point
2814 if (do_newton_visc) then
2815 strx_n = cs%newton_str_ux(i,j,qpv)
2816 stry_n = cs%newton_str_vy(i,j,qpv)
2817 strsh_n = cs%newton_str_sh(i,j,qpv)
2818 dstrain_n = (((2.*strx_n + stry_n)*ux) + ((2.*stry_n + strx_n)*vy)) + &
2819 strsh_n * (uy + vx) * 0.5
2820 endif
2821
2822 ! Newton correction for basal drag: compute inner_dot_n once per quadrature point
2823 if (use_newton) then
2824 inner_dot_n = cs%newton_umid(i,j)*uq + cs%newton_vmid(i,j)*vq
2825 endif
2826
2827 do jphi=1,2 ; jtgt = j-2+jphi ; do iphi=1,2 ; itgt = i-2+iphi
2828 if (umask(itgt,jtgt) == 1) uret_qp(iphi,jphi,qp) = ice_visc(i,j,qpv) * &
2829 (((4*ux+2*vy) * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j)) + &
2830 ((uy+vx) * phi(2*(2*(jphi-1)+iphi),qp,i,j)))
2831 if (vmask(itgt,jtgt) == 1) vret_qp(iphi,jphi,qp) = ice_visc(i,j,qpv) * &
2832 (((uy+vx) * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j)) + &
2833 ((4*vy+2*ux) * phi(2*(2*(jphi-1)+iphi),qp,i,j)))
2834
2835 ! Newton tangent stiffness correction: add (dη/dε_e^2) * (g·δε) * (g·φ_m) term
2836 if (do_newton_visc) then
2837 if (umask(itgt,jtgt) == 1) uret_qp(iphi,jphi,qp) = uret_qp(iphi,jphi,qp) + &
2838 cs%newton_visc_factor(i,j,qpv) * dstrain_n * &
2839 ((2.*strx_n + stry_n) * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j) + &
2840 strsh_n * 0.5 * phi(2*(2*(jphi-1)+iphi),qp,i,j))
2841 if (vmask(itgt,jtgt) == 1) vret_qp(iphi,jphi,qp) = vret_qp(iphi,jphi,qp) + &
2842 cs%newton_visc_factor(i,j,qpv) * dstrain_n * &
2843 (strsh_n * 0.5 * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j) + &
2844 (2.*stry_n + strx_n) * phi(2*(2*(jphi-1)+iphi),qp,i,j))
2845 endif
2846
2847 if (float_cond(i,j) == 0) then
2848 ilq = 1 ; if (iq == iphi) ilq = 2
2849 jlq = 1 ; if (jq == jphi) jlq = 2
2850 if (umask(itgt,jtgt) == 1) uret_qp(iphi,jphi,qp) = uret_qp(iphi,jphi,qp) + &
2851 ((basal_trac(i,j) * uq) * (xquad(ilq) * xquad(jlq)))
2852 if (vmask(itgt,jtgt) == 1) vret_qp(iphi,jphi,qp) = vret_qp(iphi,jphi,qp) + &
2853 ((basal_trac(i,j) * vq) * (xquad(ilq) * xquad(jlq)))
2854 ! Newton basal drag tangent stiffness: (m-1)*basal_trac/|u|^2 * u_i * (u . delta_u) contribution
2855 if (use_newton) then
2856 if (umask(itgt,jtgt) == 1) uret_qp(iphi,jphi,qp) = uret_qp(iphi,jphi,qp) + &
2857 cs%newton_drag_coef(i,j) * cs%newton_umid(i,j) * inner_dot_n * (xquad(ilq) * xquad(jlq))
2858 if (vmask(itgt,jtgt) == 1) vret_qp(iphi,jphi,qp) = vret_qp(iphi,jphi,qp) + &
2859 cs%newton_drag_coef(i,j) * cs%newton_vmid(i,j) * inner_dot_n * (xquad(ilq) * xquad(jlq))
2860 endif
2861 endif
2862 enddo ; enddo
2863 enddo ; enddo
2864
2865 !element contribution to SW node (node 1, which sees the current element as element 4)
2866 uret_b(i-1,j-1,4) = 0.25*((uret_qp(1,1,1)+uret_qp(1,1,4))+(uret_qp(1,1,2)+uret_qp(1,1,3)))
2867 vret_b(i-1,j-1,4) = 0.25*((vret_qp(1,1,1)+vret_qp(1,1,4))+(vret_qp(1,1,2)+vret_qp(1,1,3)))
2868
2869 !element contribution to NW node (node 3, which sees the current element as element 2)
2870 uret_b(i-1,j ,2) = 0.25*((uret_qp(1,2,1)+uret_qp(1,2,4))+(uret_qp(1,2,2)+uret_qp(1,2,3)))
2871 vret_b(i-1,j ,2) = 0.25*((vret_qp(1,2,1)+vret_qp(1,2,4))+(vret_qp(1,2,2)+vret_qp(1,2,3)))
2872
2873 !element contribution to SE node (node 2, which sees the current element as element 3)
2874 uret_b(i ,j-1,3) = 0.25*((uret_qp(2,1,1)+uret_qp(2,1,4))+(uret_qp(2,1,2)+uret_qp(2,1,3)))
2875 vret_b(i ,j-1,3) = 0.25*((vret_qp(2,1,1)+vret_qp(2,1,4))+(vret_qp(2,1,2)+vret_qp(2,1,3)))
2876
2877 !element contribution to NE node (node 4, which sees the current element as element 1)
2878 uret_b(i ,j ,1) = 0.25*((uret_qp(2,2,1)+uret_qp(2,2,4))+(uret_qp(2,2,2)+uret_qp(2,2,3)))
2879 vret_b(i ,j ,1) = 0.25*((vret_qp(2,2,1)+vret_qp(2,2,4))+(vret_qp(2,2,2)+vret_qp(2,2,3)))
2880
2881 if (float_cond(i,j) == 1) then
2882 ucell(:,:) = u_shlf(i-1:i,j-1:j) ; vcell(:,:) = v_shlf(i-1:i,j-1:j)
2883 hcell(:,:) = h_node(i-1:i,j-1:j)
2884
2885 call cg_action_subgrid_basal(phisub, hcell, ucell, vcell, &
2886 bathyt(i,j), dens_ratio, usub, vsub)
2887
2888 if (umask(i-1,j-1) == 1) uret_b(i-1,j-1,4) = uret_b(i-1,j-1,4) + (usub(1,1) * basal_trac(i,j))
2889 if (umask(i-1,j ) == 1) uret_b(i-1,j ,2) = uret_b(i-1,j ,2) + (usub(1,2) * basal_trac(i,j))
2890 if (umask(i ,j-1) == 1) uret_b(i ,j-1,3) = uret_b(i ,j-1,3) + (usub(2,1) * basal_trac(i,j))
2891 if (umask(i ,j ) == 1) uret_b(i ,j ,1) = uret_b(i ,j ,1) + (usub(2,2) * basal_trac(i,j))
2892
2893 if (vmask(i-1,j-1) == 1) vret_b(i-1,j-1,4) = vret_b(i-1,j-1,4) + (vsub(1,1) * basal_trac(i,j))
2894 if (vmask(i-1,j ) == 1) vret_b(i-1,j ,2) = vret_b(i-1,j ,2) + (vsub(1,2) * basal_trac(i,j))
2895 if (vmask(i ,j-1) == 1) vret_b(i ,j-1,3) = vret_b(i ,j-1,3) + (vsub(2,1) * basal_trac(i,j))
2896 if (vmask(i ,j ) == 1) vret_b(i ,j ,1) = vret_b(i ,j ,1) + (vsub(2,2) * basal_trac(i,j))
2897
2898 ! Newton basal drag correction for subgrid grounding line cells.
2899 ! inner_dot_sub(m,n) = sum over grounded sub-QPs of (u^k . delta_u) * phi_{m,n} * weight
2900 ! = newton_umid * Usub(m,n) + newton_vmid * Vsub(m,n)
2901 ! Correction to u-node (m,n): newton_drag_coef * newton_umid * inner_dot_sub(m,n)
2902 ! Correction to v-node (m,n): newton_drag_coef * newton_vmid * inner_dot_sub(m,n)
2903 if (use_newton) then
2904 if (umask(i-1,j-1)==1) uret_b(i-1,j-1,4) = uret_b(i-1,j-1,4) + cs%newton_drag_coef(i,j) * &
2905 cs%newton_umid(i,j) * ((cs%newton_umid(i,j)*usub(1,1)) + (cs%newton_vmid(i,j)*vsub(1,1)))
2906 if (umask(i-1,j )==1) uret_b(i-1,j ,2) = uret_b(i-1,j ,2) + cs%newton_drag_coef(i,j) * &
2907 cs%newton_umid(i,j) * ((cs%newton_umid(i,j)*usub(1,2)) + (cs%newton_vmid(i,j)*vsub(1,2)))
2908 if (umask(i ,j-1)==1) uret_b(i ,j-1,3) = uret_b(i ,j-1,3) + cs%newton_drag_coef(i,j) * &
2909 cs%newton_umid(i,j) * ((cs%newton_umid(i,j)*usub(2,1)) + (cs%newton_vmid(i,j)*vsub(2,1)))
2910 if (umask(i ,j )==1) uret_b(i ,j ,1) = uret_b(i ,j ,1) + cs%newton_drag_coef(i,j) * &
2911 cs%newton_umid(i,j) * ((cs%newton_umid(i,j)*usub(2,2)) + (cs%newton_vmid(i,j)*vsub(2,2)))
2912 if (vmask(i-1,j-1)==1) vret_b(i-1,j-1,4) = vret_b(i-1,j-1,4) + cs%newton_drag_coef(i,j) * &
2913 cs%newton_vmid(i,j) * ((cs%newton_umid(i,j)*usub(1,1)) + (cs%newton_vmid(i,j)*vsub(1,1)))
2914 if (vmask(i-1,j )==1) vret_b(i-1,j ,2) = vret_b(i-1,j ,2) + cs%newton_drag_coef(i,j) * &
2915 cs%newton_vmid(i,j) * ((cs%newton_umid(i,j)*usub(1,2)) + (cs%newton_vmid(i,j)*vsub(1,2)))
2916 if (vmask(i ,j-1)==1) vret_b(i ,j-1,3) = vret_b(i ,j-1,3) + cs%newton_drag_coef(i,j) * &
2917 cs%newton_vmid(i,j) * ((cs%newton_umid(i,j)*usub(2,1)) + (cs%newton_vmid(i,j)*vsub(2,1)))
2918 if (vmask(i ,j )==1) vret_b(i ,j ,1) = vret_b(i ,j ,1) + cs%newton_drag_coef(i,j) * &
2919 cs%newton_vmid(i,j) * ((cs%newton_umid(i,j)*usub(2,2)) + (cs%newton_vmid(i,j)*vsub(2,2)))
2920 endif
2921 endif
2922 endif ; enddo ; enddo
2923
2924 do j=js-1,je ; do i=is-1,ie
2925 uret(i,j) = (uret_b(i,j,1)+uret_b(i,j,4)) + (uret_b(i,j,2)+uret_b(i,j,3))
2926 vret(i,j) = (vret_b(i,j,1)+vret_b(i,j,4)) + (vret_b(i,j,2)+vret_b(i,j,3))
2927 enddo ; enddo
2928
2929end subroutine cg_action
2930
2931subroutine cg_action_subgrid_basal(Phisub, H, U, V, bathyT, dens_ratio, Ucontr, Vcontr)
2932 real, dimension(:,:,:,:,:,:), &
2933 intent(in) :: Phisub !< Quadrature structure weights at subgridscale
2934 !! locations for finite element calculations [nondim]
2935 real, dimension(2,2), intent(in) :: H !< The ice shelf thickness at nodal (corner) points [Z ~> m].
2936 real, dimension(2,2), intent(in) :: U !< The zonal ice shelf velocity at vertices [L T-1 ~> m s-1]
2937 real, dimension(2,2), intent(in) :: V !< The meridional ice shelf velocity at vertices [L T-1 ~> m s-1]
2938 real, intent(in) :: bathyT !< The depth of ocean bathymetry at tracer points
2939 !! relative to sea-level [Z ~> m].
2940 real, intent(in) :: dens_ratio !< The density of ice divided by the density
2941 !! of seawater [nondim]
2942 real, dimension(2,2), intent(out) :: Ucontr !< The areal average of u-velocities where the ice shelf
2943 !! is grounded, or 0 where it is floating [L T-1 ~> m s-1].
2944 real, dimension(2,2), intent(out) :: Vcontr !< The areal average of v-velocities where the ice shelf
2945 !! is grounded, or 0 where it is floating [L T-1 ~> m s-1].
2946
2947 real, dimension(SIZE(Phisub,3),SIZE(Phisub,3),2,2) :: Ucontr_sub, Vcontr_sub ! The contributions to Ucontr and Vcontr
2948 !! at each sub-cell
2949 real, dimension(2,2,SIZE(Phisub,3),SIZE(Phisub,3)) :: uloc_arr !The local sub-cell u-velocity [L T-1 ~> m s-1]
2950 real, dimension(2,2,SIZE(Phisub,3),SIZE(Phisub,3)) :: vloc_arr !The local sub-cell v-velocity [L T-1 ~> m s-1]
2951 real, dimension(2,2) :: Ucontr_q, Vcontr_q !Contributions to a node from each quadrature point in a sub-grid cell
2952 real :: subarea ! The fractional sub-cell area [nondim]
2953 real :: hloc ! The local sub-cell ice thickness [Z ~> m]
2954 integer :: nsub, i, j, qx, qy, m, n
2955
2956 nsub = size(phisub,3)
2957 subarea = 1.0 / (nsub**2)
2958
2959 uloc_arr(:,:,:,:) = 0.0 ; vloc_arr(:,:,:,:)=0.0
2960
2961 do j=1,nsub ; do i=1,nsub ; do qy=1,2 ; do qx=1,2
2962 hloc = ((phisub(qx,qy,i,j,1,1)*h(1,1)) + (phisub(qx,qy,i,j,2,2)*h(2,2))) + &
2963 ((phisub(qx,qy,i,j,1,2)*h(1,2)) + (phisub(qx,qy,i,j,2,1)*h(2,1)))
2964 if (dens_ratio * hloc - bathyt > 0) then
2965 uloc_arr(qx,qy,i,j) = (((phisub(qx,qy,i,j,1,1) * u(1,1)) + (phisub(qx,qy,i,j,2,2) * u(2,2))) + &
2966 ((phisub(qx,qy,i,j,1,2) * u(1,2)) + (phisub(qx,qy,i,j,2,1) * u(2,1))))
2967 vloc_arr(qx,qy,i,j) = (((phisub(qx,qy,i,j,1,1) * v(1,1)) + (phisub(qx,qy,i,j,2,2) * v(2,2))) + &
2968 ((phisub(qx,qy,i,j,1,2) * v(1,2)) + (phisub(qx,qy,i,j,2,1) * v(2,1))))
2969 endif
2970 enddo ; enddo ; enddo ; enddo
2971
2972 do n=1,2 ; do m=1,2 ; do j=1,nsub ; do i=1,nsub
2973 do qy=1,2 ; do qx=1,2
2974 !calculate quadrature point contributions for the sub-cell, to each node
2975 ucontr_q(qx,qy) = phisub(qx,qy,i,j,m,n) * uloc_arr(qx,qy,i,j)
2976 vcontr_q(qx,qy) = phisub(qx,qy,i,j,m,n) * vloc_arr(qx,qy,i,j)
2977 enddo ; enddo
2978
2979 !calculate sub-cell contribution to each node by summing up quadrature point contributions from the sub-cell
2980 ucontr_sub(i,j,m,n) = (subarea * 0.25) * ((ucontr_q(1,1) + ucontr_q(2,2)) + (ucontr_q(1,2)+ucontr_q(2,1)))
2981 vcontr_sub(i,j,m,n) = (subarea * 0.25) * ((vcontr_q(1,1) + vcontr_q(2,2)) + (vcontr_q(1,2)+vcontr_q(2,1)))
2982 enddo ; enddo ; enddo ; enddo
2983
2984 !sum up the sub-cell contributions to each node
2985 do n=1,2 ; do m=1,2
2986 call sum_square_matrix(ucontr(m,n),ucontr_sub(:,:,m,n),nsub)
2987 call sum_square_matrix(vcontr(m,n),vcontr_sub(:,:,m,n),nsub)
2988 enddo ; enddo
2989
2990end subroutine cg_action_subgrid_basal
2991
2992
2993!! Returns the sum of the elements in a square matrix. This sum is bitwise identical even if the matrices are rotated.
2994subroutine sum_square_matrix(sum_out, mat_in, n)
2995 integer, intent(in) :: n !< The length and width of each matrix in mat_in
2996 real, dimension(n,n), intent(in) :: mat_in !< The n x n matrix whose elements will be summed
2997 real, intent(out) :: sum_out !< The sum of the elements of matrix mat_in
2998 integer :: s0, e0, s1, e1
2999
3000 sum_out = 0.0
3001
3002 s0 = 1 ; e0 = n
3003
3004 !start by summing elements on outer edges of matrix
3005 do while (s0<e0)
3006
3007 !corners
3008 sum_out = sum_out + ( (mat_in(s0,s0) + mat_in(e0,e0)) + (mat_in(e0,s0) + mat_in(s0,e0)) )
3009
3010 s1 = s0+1 ; e1 = e0-1
3011
3012 do while (s1<e1) !non-corners
3013
3014 sum_out = sum_out + &
3015 ( ( (mat_in(s0,s1) + mat_in(s1,s0)) + (mat_in(e0,e1) + mat_in(e1,e0)) ) + &
3016 ( (mat_in(e1,s0) + mat_in(e0,s1)) + (mat_in(s1,e0) + mat_in(s0,e1)) ) )
3017
3018 s1 = s1+1 ; e1 = e1-1
3019 enddo
3020
3021 !center element of an edge
3022 if (s1==e1) sum_out = sum_out + ( (mat_in(s1,s0) + mat_in(e1,e0)) + (mat_in(e0,e1) + mat_in(s0,s1)) )
3023
3024 s0 = s0+1 ; e0 = e0-1 !next loop iteration using new edges that are one element inward of the current edges
3025 enddo
3026
3027 !center element of entire matrix
3028 if (s0==e0) sum_out = sum_out + mat_in(s0,e0)
3029
3030end subroutine sum_square_matrix
3031
3032!> returns the diagonal entries of the matrix for a Jacobi preconditioning
3033subroutine matrix_diagonal(CS, G, US, float_cond, H_node, ice_visc, basal_trac, hmask, dens_ratio, &
3034 Phi, Phisub, u_diagonal, v_diagonal)
3035
3036 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
3037 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
3038 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
3039 real, dimension(SZDI_(G),SZDJ_(G)), &
3040 intent(in) :: float_cond !< If GL_regularize=true, indicates cells containing
3041 !! the grounding line (float_cond=1) or not (float_cond=0)
3042 real, dimension(SZDIB_(G),SZDJB_(G)), &
3043 intent(in) :: H_node !< The ice shelf thickness at nodal
3044 !! (corner) points [Z ~> m].
3045 real, dimension(SZDI_(G),SZDJ_(G),CS%visc_qps), &
3046 intent(in) :: ice_visc !< A field related to the ice viscosity from Glen's
3047 !! flow law [R L4 Z T-1 ~> kg m2 s-1].
3048 real, dimension(SZDI_(G),SZDJ_(G)), &
3049 intent(in) :: basal_trac !< Area-integrated taub_beta field related to the nonlinear
3050 !! part of the "linearized" basal stress [R L3 T-1 ~> kg s-1].
3051 real, dimension(SZDI_(G),SZDJ_(G)), &
3052 intent(in) :: hmask !< A mask indicating which tracer points are
3053 !! partly or fully covered by an ice-shelf
3054 real, intent(in) :: dens_ratio !< The density of ice divided by the density
3055 !! of seawater [nondim]
3056 real, dimension(8,4,SZDI_(G),SZDJ_(G)), &
3057 intent(in) :: Phi !< The gradients of bilinear basis elements at Gaussian
3058 !! quadrature points surrounding the cell vertices [L-1 ~> m-1]
3059 real, dimension(:,:,:,:,:,:), intent(in) :: Phisub !< Quadrature structure weights at subgridscale
3060 !! locations for finite element calculations [nondim]
3061 real, dimension(SZDIB_(G),SZDJB_(G)), &
3062 intent(inout) :: u_diagonal !< The diagonal elements of the u-velocity
3063 !! matrix from the left-hand side of the solver [R L2 Z T-1 ~> kg s-1]
3064 real, dimension(SZDIB_(G),SZDJB_(G)), &
3065 intent(inout) :: v_diagonal !< The diagonal elements of the v-velocity
3066 !! matrix from the left-hand side of the solver [R L2 Z T-1 ~> kg s-1]
3067
3068
3069! returns the diagonal entries of the matrix for a Jacobi preconditioning
3070
3071 real :: ux, uy, vx, vy ! Interpolated weight gradients [L-1 ~> m-1]
3072 real :: uq, vq
3073 real :: strx_n, stry_n, strsh_n ! Newton viscosity strain rates [T-1 ~> s-1]
3074 real :: dstrain_diag_u, dstrain_diag_v ! Newton viscosity diagonal correction factors [T-1 L-1 ~> s-1 m-1]
3075 real :: phi_m_sq ! Squared basis function value at quadrature point [nondim]
3076 real, dimension(2) :: xquad
3077 real, dimension(2,2) :: Hcell, sub_ground
3078 real, dimension(2,2,4) :: u_diag_qp, v_diag_qp
3079 real, dimension(SZDIB_(G),SZDJB_(G),4) :: u_diag_b, v_diag_b
3080 logical :: do_newton_visc ! Whether to apply viscosity-related Newton tangent stiffness corrections
3081 logical :: visc_qp4
3082 integer :: i, j, isc, jsc, iec, jec, iphi, jphi, iq, jq, ilq, jlq, Itgt, Jtgt, qp, qpv
3083
3084 isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
3085
3086 xquad(1) = .5 * (1-sqrt(1./3)) ; xquad(2) = .5 * (1+sqrt(1./3))
3087
3088 if (cs%visc_qps == 4) then
3089 visc_qp4=.true.
3090 else
3091 visc_qp4=.false.
3092 qpv = 1
3093 endif
3094
3095 do_newton_visc = cs%doing_newton .and. trim(cs%ice_viscosity_compute) == "MODEL"
3096
3097 u_diag_b(:,:,:)=0.0
3098 v_diag_b(:,:,:)=0.0
3099
3100 do j=jsc-1,jec+1 ; do i=isc-1,iec+1 ; if (hmask(i,j) == 1 .or. hmask(i,j)==3) then
3101
3102 ! Phi(2*i-1,j) gives d(Phi_i)/dx at quadrature point j
3103 ! Phi(2*i,j) gives d(Phi_i)/dy at quadrature point j
3104
3105 u_diag_qp(:,:,:) = 0.0 ; v_diag_qp(:,:,:) = 0.0
3106
3107 do iq=1,2 ; do jq=1,2
3108
3109 qp = 2*(jq-1)+iq !current quad point
3110 if (visc_qp4) qpv = qp !current quad point for viscosity
3111
3112 ! Pre-compute Newton strain data for this QP (for viscosity diagonal correction)
3113 if (do_newton_visc) then
3114 strx_n = cs%newton_str_ux(i,j,qpv)
3115 stry_n = cs%newton_str_vy(i,j,qpv)
3116 strsh_n = cs%newton_str_sh(i,j,qpv)
3117 endif
3118
3119 do jphi=1,2 ; jtgt = j-2+jphi ; do iphi=1,2 ; itgt = i-2+iphi
3120
3121 ilq = 1 ; if (iq == iphi) ilq = 2
3122 jlq = 1 ; if (jq == jphi) jlq = 2
3123 phi_m_sq = (xquad(ilq) * xquad(jlq))**2
3124
3125 if (cs%umask(itgt,jtgt) == 1) then
3126
3127 ux = phi(2*(2*(jphi-1)+iphi)-1,qp,i,j)
3128 uy = phi(2*(2*(jphi-1)+iphi),qp,i,j)
3129 vx = 0.
3130 vy = 0.
3131
3132 u_diag_qp(iphi,jphi,qp) = &
3133 ice_visc(i,j,qpv) * (((4*ux+2*vy) * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j)) + &
3134 ((uy+vx) * phi(2*(2*(jphi-1)+iphi),qp,i,j)))
3135
3136 ! Newton viscosity diagonal correction: newton_visc_factor * (g . grad_phi_m_u)^2
3137 ! where grad_phi_m_u = [(2*strx+stry)*Phi_xm + strsh/2*Phi_ym] for u-DOF at node m
3138 if (do_newton_visc) then
3139 dstrain_diag_u = (2.*strx_n + stry_n) * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j) + &
3140 strsh_n * 0.5 * phi(2*(2*(jphi-1)+iphi),qp,i,j)
3141 u_diag_qp(iphi,jphi,qp) = u_diag_qp(iphi,jphi,qp) + &
3142 cs%newton_visc_factor(i,j,qpv) * dstrain_diag_u**2
3143 endif
3144
3145 if (float_cond(i,j) == 0) then
3146 uq = xquad(ilq) * xquad(jlq)
3147 u_diag_qp(iphi,jphi,qp) = u_diag_qp(iphi,jphi,qp) + &
3148 (basal_trac(i,j) * uq) * (xquad(ilq) * xquad(jlq))
3149 ! Newton basal drag diagonal correction: newton_drag_coef * (umid_i)^2 * phi_m^2
3150 if (cs%doing_newton) then
3151 u_diag_qp(iphi,jphi,qp) = u_diag_qp(iphi,jphi,qp) + &
3152 cs%newton_drag_coef(i,j) * cs%newton_umid(i,j)**2 * phi_m_sq
3153 endif
3154 endif
3155 endif
3156
3157 if (cs%vmask(itgt,jtgt) == 1) then
3158
3159 vx = phi(2*(2*(jphi-1)+iphi)-1,qp,i,j)
3160 vy = phi(2*(2*(jphi-1)+iphi),qp,i,j)
3161 ux = 0.
3162 uy = 0.
3163
3164 v_diag_qp(iphi,jphi,qp) = &
3165 ice_visc(i,j,qpv) * (((uy+vx) * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j)) + &
3166 ((4*vy+2*ux) * phi(2*(2*(jphi-1)+iphi),qp,i,j)))
3167
3168 ! Newton viscosity diagonal correction for v-DOF: uses [strsh/2*Phi_xm + (2*stry+strx)*Phi_ym]
3169 if (do_newton_visc) then
3170 dstrain_diag_v = strsh_n * 0.5 * phi(2*(2*(jphi-1)+iphi)-1,qp,i,j) + &
3171 (2.*stry_n + strx_n) * phi(2*(2*(jphi-1)+iphi),qp,i,j)
3172 v_diag_qp(iphi,jphi,qp) = v_diag_qp(iphi,jphi,qp) + &
3173 cs%newton_visc_factor(i,j,qpv) * dstrain_diag_v**2
3174 endif
3175
3176 if (float_cond(i,j) == 0) then
3177 vq = xquad(ilq) * xquad(jlq)
3178 v_diag_qp(iphi,jphi,qp) = v_diag_qp(iphi,jphi,qp) + &
3179 (basal_trac(i,j) * vq) * (xquad(ilq) * xquad(jlq))
3180 ! Newton basal drag diagonal correction: newton_drag_coef * (vmid_i)^2 * phi_m^2
3181 if (cs%doing_newton) then
3182 v_diag_qp(iphi,jphi,qp) = v_diag_qp(iphi,jphi,qp) + &
3183 cs%newton_drag_coef(i,j) * cs%newton_vmid(i,j)**2 * phi_m_sq
3184 endif
3185 endif
3186 endif
3187 enddo ; enddo
3188 enddo ; enddo
3189
3190 !element contribution to SW node (node 1, which sees the current element as element 4)
3191 u_diag_b(i-1,j-1,4) = 0.25*((u_diag_qp(1,1,1)+u_diag_qp(1,1,4))+(u_diag_qp(1,1,2)+u_diag_qp(1,1,3)))
3192 v_diag_b(i-1,j-1,4) = 0.25*((v_diag_qp(1,1,1)+v_diag_qp(1,1,4))+(v_diag_qp(1,1,2)+v_diag_qp(1,1,3)))
3193
3194 !element contribution to NW node (node 3, which sees the current element as element 2)
3195 u_diag_b(i-1,j ,2) = 0.25*((u_diag_qp(1,2,1)+u_diag_qp(1,2,4))+(u_diag_qp(1,2,2)+u_diag_qp(1,2,3)))
3196 v_diag_b(i-1,j ,2) = 0.25*((v_diag_qp(1,2,1)+v_diag_qp(1,2,4))+(v_diag_qp(1,2,2)+v_diag_qp(1,2,3)))
3197
3198 !element contribution to SE node (node 2, which sees the current element as element 3)
3199 u_diag_b(i ,j-1,3) = 0.25*((u_diag_qp(2,1,1)+u_diag_qp(2,1,4))+(u_diag_qp(2,1,2)+u_diag_qp(2,1,3)))
3200 v_diag_b(i ,j-1,3) = 0.25*((v_diag_qp(2,1,1)+v_diag_qp(2,1,4))+(v_diag_qp(2,1,2)+v_diag_qp(2,1,3)))
3201
3202 !element contribution to NE node (node 4, which sees the current element as element 1)
3203 u_diag_b(i ,j ,1) = 0.25*((u_diag_qp(2,2,1)+u_diag_qp(2,2,4))+(u_diag_qp(2,2,2)+u_diag_qp(2,2,3)))
3204 v_diag_b(i ,j ,1) = 0.25*((v_diag_qp(2,2,1)+v_diag_qp(2,2,4))+(v_diag_qp(2,2,2)+v_diag_qp(2,2,3)))
3205
3206 if (float_cond(i,j) == 1) then
3207 hcell(:,:) = h_node(i-1:i,j-1:j)
3208 call cg_diagonal_subgrid_basal(phisub, hcell, cs%bed_elev(i,j), dens_ratio, sub_ground)
3209
3210 if (cs%umask(i-1,j-1) == 1) u_diag_b(i-1,j-1,4) = u_diag_b(i-1,j-1,4) + (sub_ground(1,1) * basal_trac(i,j))
3211 if (cs%umask(i-1,j ) == 1) u_diag_b(i-1,j ,2) = u_diag_b(i-1,j ,2) + (sub_ground(1,2) * basal_trac(i,j))
3212 if (cs%umask(i ,j-1) == 1) u_diag_b(i ,j-1,3) = u_diag_b(i ,j-1,3) + (sub_ground(2,1) * basal_trac(i,j))
3213 if (cs%umask(i ,j ) == 1) u_diag_b(i ,j ,1) = u_diag_b(i ,j ,1) + (sub_ground(2,2) * basal_trac(i,j))
3214
3215 if (cs%vmask(i-1,j-1) == 1) v_diag_b(i-1,j-1,4) = v_diag_b(i-1,j-1,4) + (sub_ground(1,1) * basal_trac(i,j))
3216 if (cs%vmask(i-1,j ) == 1) v_diag_b(i-1,j ,2) = v_diag_b(i-1,j ,2) + (sub_ground(1,2) * basal_trac(i,j))
3217 if (cs%vmask(i ,j-1) == 1) v_diag_b(i ,j-1,3) = v_diag_b(i ,j-1,3) + (sub_ground(2,1) * basal_trac(i,j))
3218 if (cs%vmask(i ,j ) == 1) v_diag_b(i ,j ,1) = v_diag_b(i ,j ,1) + (sub_ground(2,2) * basal_trac(i,j))
3219
3220 ! Newton basal drag diagonal correction for subgrid grounding line cells.
3221 ! sub_ground(m,n) = sum over grounded sub-QPs of phi_{m,n}^2 * weight, computed by
3222 ! CG_diagonal_subgrid_basal. Newton diagonal = newton_drag_coef * umid^2 * sub_ground (for u-block).
3223 if (cs%doing_newton) then
3224 if (cs%umask(i-1,j-1)==1) u_diag_b(i-1,j-1,4) = u_diag_b(i-1,j-1,4) + &
3225 cs%newton_drag_coef(i,j) * cs%newton_umid(i,j)**2 * sub_ground(1,1)
3226 if (cs%umask(i-1,j )==1) u_diag_b(i-1,j ,2) = u_diag_b(i-1,j ,2) + &
3227 cs%newton_drag_coef(i,j) * cs%newton_umid(i,j)**2 * sub_ground(1,2)
3228 if (cs%umask(i ,j-1)==1) u_diag_b(i ,j-1,3) = u_diag_b(i ,j-1,3) + &
3229 cs%newton_drag_coef(i,j) * cs%newton_umid(i,j)**2 * sub_ground(2,1)
3230 if (cs%umask(i ,j )==1) u_diag_b(i ,j ,1) = u_diag_b(i ,j ,1) + &
3231 cs%newton_drag_coef(i,j) * cs%newton_umid(i,j)**2 * sub_ground(2,2)
3232 if (cs%vmask(i-1,j-1)==1) v_diag_b(i-1,j-1,4) = v_diag_b(i-1,j-1,4) + &
3233 cs%newton_drag_coef(i,j) * cs%newton_vmid(i,j)**2 * sub_ground(1,1)
3234 if (cs%vmask(i-1,j )==1) v_diag_b(i-1,j ,2) = v_diag_b(i-1,j ,2) + &
3235 cs%newton_drag_coef(i,j) * cs%newton_vmid(i,j)**2 * sub_ground(1,2)
3236 if (cs%vmask(i ,j-1)==1) v_diag_b(i ,j-1,3) = v_diag_b(i ,j-1,3) + &
3237 cs%newton_drag_coef(i,j) * cs%newton_vmid(i,j)**2 * sub_ground(2,1)
3238 if (cs%vmask(i ,j )==1) v_diag_b(i ,j ,1) = v_diag_b(i ,j ,1) + &
3239 cs%newton_drag_coef(i,j) * cs%newton_vmid(i,j)**2 * sub_ground(2,2)
3240 endif
3241 endif
3242 endif ; enddo ; enddo
3243
3244 do j=jsc-2,jec+1 ; do i=isc-2,iec+1
3245 u_diagonal(i,j) = (u_diag_b(i,j,1)+u_diag_b(i,j,4)) + (u_diag_b(i,j,2)+u_diag_b(i,j,3))
3246 v_diagonal(i,j) = (v_diag_b(i,j,1)+v_diag_b(i,j,4)) + (v_diag_b(i,j,2)+v_diag_b(i,j,3))
3247 enddo ; enddo
3248
3249end subroutine matrix_diagonal
3250
3251subroutine cg_diagonal_subgrid_basal (Phisub, H_node, bathyT, dens_ratio, f_grnd)
3252 real, dimension(:,:,:,:,:,:), &
3253 intent(in) :: Phisub !< Quadrature structure weights at subgridscale
3254 !! locations for finite element calculations [nondim]
3255 real, dimension(2,2), intent(in) :: H_node !< The ice shelf thickness at nodal (corner)
3256 !! points [Z ~> m].
3257 real, intent(in) :: bathyT !< The depth of ocean bathymetry at tracer points [Z ~> m].
3258 real, intent(in) :: dens_ratio !< The density of ice divided by the density
3259 !! of seawater [nondim]
3260 real, dimension(2,2), intent(out) :: f_grnd !< The weighted fraction of the sub-cell where the ice shelf
3261 !! is grounded [nondim]
3262
3263 real, dimension(SIZE(Phisub,3),SIZE(Phisub,3),2,2) :: f_grnd_sub ! The contributions to nodal f_grnd
3264 !! from each sub-cell
3265 integer, dimension(2,2,SIZE(Phisub,3),SIZE(Phisub,3)) :: grnd_stat !0 at floating quad points, 1 at grounded
3266 real, dimension(2,2) :: f_grnd_q !Contributions to a node from each quadrature point in a sub-grid cell
3267 real :: subarea ! The fractional sub-cell area [nondim]
3268 real :: hloc ! The local sub-region thickness [Z ~> m]
3269 integer :: nsub, i, j, qx, qy, m, n
3270
3271 nsub = size(phisub,3)
3272 subarea = 1.0 / (nsub**2)
3273
3274 grnd_stat(:,:,:,:)=0
3275
3276 do j=1,nsub ; do i=1,nsub ; do qy=1,2 ; do qx=1,2
3277 hloc = ((phisub(qx,qy,i,j,1,1)*h_node(1,1)) + (phisub(qx,qy,i,j,2,2)*h_node(2,2))) + &
3278 ((phisub(qx,qy,i,j,1,2)*h_node(1,2)) + (phisub(qx,qy,i,j,2,1)*h_node(2,1)))
3279 if (dens_ratio * hloc - bathyt > 0) grnd_stat(qx,qy,i,j) = 1
3280 enddo ; enddo ; enddo ; enddo
3281
3282 do n=1,2 ; do m=1,2 ; do j=1,nsub ; do i=1,nsub
3283 do qy=1,2 ; do qx = 1,2
3284 f_grnd_q(qx,qy) = grnd_stat(qx,qy,i,j) * phisub(qx,qy,i,j,m,n)**2
3285 enddo ; enddo
3286 !calculate sub-cell contribution to each node by summing up quadrature point contributions from the sub-cell
3287 f_grnd_sub(i,j,m,n) = (subarea * 0.25) * ((f_grnd_q(1,1) + f_grnd_q(2,2)) + (f_grnd_q(1,2)+f_grnd_q(2,1)))
3288 enddo ; enddo ; enddo ; enddo
3289
3290 !sum up the sub-cell contributions to each node
3291 do n=1,2 ; do m=1,2
3292 call sum_square_matrix(f_grnd(m,n),f_grnd_sub(:,:,m,n),nsub)
3293 enddo ; enddo
3294
3295end subroutine cg_diagonal_subgrid_basal
3296
3297!> Post_data calls related to ice-sheet flux divergence, strain-rate, and deviatoric stress
3298subroutine is_dynamics_post_data_2(CS, ISS, G)
3299 type(ice_shelf_dyn_cs), intent(inout) :: CS !< A pointer to the ice shelf control structure
3300 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
3301 !! the ice-shelf state
3302 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
3303 real, dimension(SZDIB_(G),SZDJB_(G)) :: H_node ! Ice shelf thickness at corners [Z ~> m].
3304 real, dimension(SZDIB_(G),SZDJB_(G)) :: Hu ! Ice shelf u_flux at corners [Z L T-1 ~> m2 s-1].
3305 real, dimension(SZDIB_(G),SZDJB_(G)) :: Hv ! Ice shelf v_flux at corners [Z L T-1 ~> m2 s-1].
3306 real, dimension(SZDI_(G),SZDJ_(G)) :: Hux ! Ice shelf d(u_flux)/dx at cell centers [Z T-1 ~> m s-1].
3307 real, dimension(SZDI_(G),SZDJ_(G)) :: Hvy ! Ice shelf d(v_flux)/dy at cell centers [Z T-1 ~> m s-1].
3308 real, dimension(SZDI_(G),SZDJ_(G)) :: flux_div ! horizontal flux divergence div(uH) [Z T-1 ~> m s-1].
3309 real, dimension(SZDI_(G),SZDJ_(G),3) :: strain_rate ! strain-rate components xx,yy, and xy [T-1 ~> s-1]
3310 real, dimension(SZDI_(G),SZDJ_(G),2) :: p_strain_rate ! horizontal principal strain-rates [T-1 ~> s-1]
3311 real, dimension(SZDI_(G),SZDJ_(G),3) :: dev_stress ! deviatoric stress components xx,yy, and xy [R L Z T-2 ~> Pa]
3312 real, dimension(SZDI_(G),SZDJ_(G),2) :: p_dev_stress ! horizontal principal deviatoric stress [R L Z T-2 ~> Pa]
3313 real, dimension(SZDI_(G),SZDJ_(G)) :: ice_visc ! area-averaged ice viscosity [R L2 T-1 ~> Pa s]
3314 real :: p1,p2 ! Used to calculate strain-rate principal components [T-1 ~> s-1]
3315 integer :: i, j
3316
3317 !Allocate the gradient basis functions for 1 cell-centered quadrature point per cell
3318 if (.not. associated(cs%PhiC)) then
3319 allocate(cs%PhiC(1:8,g%isc:g%iec,g%jsc:g%jec), source=0.0)
3320 do j=g%jsc,g%jec ; do i=g%isc,g%iec
3321 call bilinear_shape_fn_grid_1qp(g, i, j, cs%PhiC(:,i,j))
3322 enddo ; enddo
3323 endif
3324
3325 !Calculate flux divergence and its components
3326 if (cs%id_duHdx > 0 .or. cs%id_dvHdy > 0 .or. cs%id_fluxdiv > 0) then
3327 call interpolate_h_to_b(g, iss%h_shelf, iss%hmask, h_node, cs%min_h_shelf)
3328
3329 hu(:,:) = 0.0 ; hv(:,:) = 0.0 ; hux(:,:) = 0.0 ; hvy(:,:) = 0.0 ; flux_div(:,:) = 0.0
3330 do j=g%jscB,g%jecB ; do i=g%iscB,g%iecB
3331 if (cs%umask(i,j) > 0) then
3332 hu(i,j) = (h_node(i,j) * cs%u_shelf(i,j))
3333 endif
3334 if (cs%vmask(i,j) > 0) then
3335 hv(i,j) = (h_node(i,j) * cs%v_shelf(i,j))
3336 endif
3337 enddo ; enddo
3338
3339 do j=g%jsc,g%jec ; do i=g%isc,g%iec
3340 if ((iss%hmask(i,j) == 1) .or. (iss%hmask(i,j) == 3)) then
3341 !components of flux divergence at cell centers
3342 hux(i,j) = (((hu(i-1,j-1) * cs%PhiC(1,i,j)) + (hu(i,j ) * cs%PhiC(7,i,j))) + &
3343 ((hu(i-1,j ) * cs%PhiC(5,i,j)) + (hu(i,j-1) * cs%PhiC(3,i,j))))
3344
3345 hvy(i,j) = (((hv(i-1,j-1) * cs%PhiC(2,i,j)) + (hv(i,j ) * cs%PhiC(8,i,j))) + &
3346 ((hv(i-1,j ) * cs%PhiC(6,i,j)) + (hv(i,j-1) * cs%PhiC(4,i,j))))
3347 flux_div(i,j) = hux(i,j) + hvy(i,j)
3348 endif
3349 enddo ; enddo
3350
3351 if (cs%id_duHdx > 0) call post_data(cs%id_duHdx, hux, cs%diag)
3352 if (cs%id_dvHdy > 0) call post_data(cs%id_dvHdy, hvy, cs%diag)
3353 if (cs%id_fluxdiv > 0) call post_data(cs%id_fluxdiv, flux_div, cs%diag)
3354 endif
3355
3356 if (cs%id_devstress_xx > 0 .or. cs%id_devstress_yy > 0 .or. cs%id_devstress_xy > 0 .or. &
3357 cs%id_strainrate_xx > 0 .or. cs%id_strainrate_yy > 0 .or. cs%id_strainrate_xy > 0 .or. &
3358 cs%id_pdevstress_1 > 0 .or. cs%id_pdevstress_2 > 0 .or. &
3359 cs%id_pstrainrate_1 > 0 .or. cs%id_pstrainrate_2 > 0) then
3360
3361 strain_rate(:,:,:) = 0.0
3362 do j=g%jsc,g%jec ; do i=g%isc,g%iec
3363 !strain-rates at cell centers
3364 if ((iss%hmask(i,j) == 1) .or. (iss%hmask(i,j) == 3)) then
3365 !strain_rate(:,:,1) = strain_rate_xx(:,:) = ux(:,:)
3366 strain_rate(i,j,1) = (((cs%u_shelf(i-1,j-1) * cs%PhiC(1,i,j)) + (cs%u_shelf(i,j ) * cs%PhiC(7,i,j))) + &
3367 ((cs%u_shelf(i-1,j ) * cs%PhiC(5,i,j)) + (cs%u_shelf(i,j-1) * cs%PhiC(3,i,j))))
3368 !strain_rate(:,:,2) = strain_rate_yy(:,:) = uy(:,:)
3369 strain_rate(i,j,2) = (((cs%v_shelf(i-1,j-1) * cs%PhiC(2,i,j)) + (cs%v_shelf(i,j ) * cs%PhiC(8,i,j))) + &
3370 ((cs%v_shelf(i-1,j ) * cs%PhiC(6,i,j)) + (cs%v_shelf(i,j-1) * cs%PhiC(4,i,j))))
3371 !strain_rate(:,:,3) = strain_rate_xy(:,:) = 0.5 * (uy(:,:) + vy(:,:))
3372 strain_rate(i,j,3) = 0.5 * ((((cs%u_shelf(i-1,j-1) * cs%PhiC(2,i,j)) + (cs%u_shelf(i,j ) * cs%PhiC(8,i,j))) + &
3373 ((cs%u_shelf(i-1,j ) * cs%PhiC(6,i,j)) + (cs%u_shelf(i,j-1) * cs%PhiC(4,i,j))))+ &
3374 (((cs%v_shelf(i-1,j-1) * cs%PhiC(1,i,j)) + (cs%v_shelf(i,j ) * cs%PhiC(7,i,j))) + &
3375 ((cs%v_shelf(i-1,j ) * cs%PhiC(5,i,j)) + (cs%v_shelf(i,j-1) * cs%PhiC(3,i,j)))))
3376 endif
3377 enddo ; enddo
3378
3379
3380 if (cs%id_strainrate_xx > 0) call post_data(cs%id_strainrate_xx, strain_rate(:,:,1), cs%diag)
3381 if (cs%id_strainrate_yy > 0) call post_data(cs%id_strainrate_yy, strain_rate(:,:,2), cs%diag)
3382 if (cs%id_strainrate_xy > 0) call post_data(cs%id_strainrate_xy, strain_rate(:,:,3), cs%diag)
3383
3384 if (cs%id_pstrainrate_1 > 0 .or. cs%id_pstrainrate_2 > 0 .or. &
3385 cs%id_pdevstress_1 > 0 .or. cs%id_pdevstress_2 > 0) then
3386 p_strain_rate(:,:,:) = 0.0
3387 do j=g%jsc,g%jec ; do i=g%isc,g%iec
3388 p1 = 0.5*( strain_rate(i,j,1) + strain_rate(i,j,2))
3389 p2 = sqrt( (( 0.5 * (strain_rate(i,j,1) - strain_rate(i,j,2)) )**2) + (strain_rate(i,j,3)**2) )
3390 p_strain_rate(i,j,1) = p1+p2 !Max horizontal principal strain-rate
3391 p_strain_rate(i,j,2) = p1-p2 !Min horizontal principal strain-rate
3392 enddo ; enddo
3393
3394 if (cs%id_pstrainrate_1 > 0) call post_data(cs%id_pstrainrate_1, p_strain_rate(:,:,1), cs%diag)
3395 if (cs%id_pstrainrate_2 > 0) call post_data(cs%id_pstrainrate_2, p_strain_rate(:,:,2), cs%diag)
3396 endif
3397
3398 if (cs%id_devstress_xx > 0 .or. cs%id_devstress_yy > 0 .or. cs%id_devstress_xy > 0 .or. &
3399 cs%id_pdevstress_1 > 0 .or. cs%id_pdevstress_2 > 0) then
3400
3401 call ice_visc_diag(cs,g,ice_visc)
3402
3403 if (cs%id_devstress_xx > 0 .or. cs%id_devstress_yy > 0 .or. cs%id_devstress_xy > 0) then
3404 dev_stress(:,:,:)=0.0
3405 do j=g%jsc,g%jec ; do i=g%isc,g%iec
3406 if (iss%h_shelf(i,j)>0) then
3407 dev_stress(i,j,1) = 2*ice_visc(i,j)*strain_rate(i,j,1)/iss%h_shelf(i,j) !deviatoric stress xx
3408 dev_stress(i,j,2) = 2*ice_visc(i,j)*strain_rate(i,j,2)/iss%h_shelf(i,j) !deviatoric stress yy
3409 dev_stress(i,j,3) = 2*ice_visc(i,j)*strain_rate(i,j,3)/iss%h_shelf(i,j) !deviatoric stress xy
3410 endif
3411 enddo ; enddo
3412 if (cs%id_devstress_xx > 0) call post_data(cs%id_devstress_xx, dev_stress(:,:,1), cs%diag)
3413 if (cs%id_devstress_yy > 0) call post_data(cs%id_devstress_yy, dev_stress(:,:,2), cs%diag)
3414 if (cs%id_devstress_xy > 0) call post_data(cs%id_devstress_xy, dev_stress(:,:,3), cs%diag)
3415 endif
3416
3417 if (cs%id_pdevstress_1 > 0 .or. cs%id_pdevstress_2 > 0) then
3418 p_dev_stress(:,:,:)=0.0
3419 do j=g%jsc,g%jec ; do i=g%isc,g%iec
3420 if (iss%h_shelf(i,j)>0) then
3421 p_dev_stress(i,j,1) = 2*ice_visc(i,j)*p_strain_rate(i,j,1)/iss%h_shelf(i,j) !max horiz principal dev stress
3422 p_dev_stress(i,j,2) = 2*ice_visc(i,j)*p_strain_rate(i,j,2)/iss%h_shelf(i,j) !min horiz principal dev stress
3423 endif
3424 enddo ; enddo
3425 if (cs%id_pdevstress_1 > 0) call post_data(cs%id_pdevstress_1, p_dev_stress(:,:,1), cs%diag)
3426 if (cs%id_pdevstress_2 > 0) call post_data(cs%id_pdevstress_2, p_dev_stress(:,:,2), cs%diag)
3427 endif
3428 endif
3429 endif
3430end subroutine is_dynamics_post_data_2
3431
3432!> Update depth integrated viscosity, based on horizontal strain rates
3433subroutine calc_shelf_visc(CS, ISS, G, US, u_shlf, v_shlf)
3434 type(ice_shelf_dyn_cs), intent(inout) :: CS !< A pointer to the ice shelf control structure
3435 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
3436 !! the ice-shelf state
3437 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
3438 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
3439 real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
3440 intent(inout) :: u_shlf !< The zonal ice shelf velocity [L T-1 ~> m s-1].
3441 real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
3442 intent(inout) :: v_shlf !< The meridional ice shelf velocity [L T-1 ~> m s-1].
3443
3444! update DEPTH_INTEGRATED viscosity, based on horizontal strain rates - this is for bilinear FEM solve
3445
3446
3447! this may be subject to change later... to make it "hybrid"
3448! real, dimension(SZDIB_(G),SZDJB_(G)) :: eII, ux, uy, vx, vy
3449 integer :: i, j, iscq, iecq, jscq, jecq, isd, jsd, ied, jed, iegq, jegq, iq, jq
3450 integer :: giec, gjec, gisc, gjsc, isc, jsc, iec, jec, is, js
3451 real :: Visc_coef, n_g
3452 real :: ux, uy, vx, vy
3453 real :: eps_min ! Velocity shears [T-1 ~> s-1]
3454 logical :: model_qp1, model_qp4
3455
3456 isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
3457 iscq = g%iscB ; iecq = g%iecB ; jscq = g%jscB ; jecq = g%jecB
3458 isd = g%isd ; jsd = g%jsd ; ied = g%ied ; jed = g%jed
3459 iegq = g%iegB ; jegq = g%jegB
3460 gisc = g%domain%nihalo+1 ; gjsc = g%domain%njhalo+1
3461 giec = g%domain%niglobal+gisc ; gjec = g%domain%njglobal+gjsc
3462 is = iscq - 1 ; js = jscq - 1
3463
3464 if (trim(cs%ice_viscosity_compute) == "MODEL") then
3465 if (cs%visc_qps==1) then
3466 model_qp1=.true.
3467 model_qp4=.false.
3468 else
3469 model_qp1=.false.
3470 model_qp4=.true.
3471 endif
3472 endif
3473
3474 n_g = cs%n_glen ; eps_min = cs%eps_glen_min
3475
3476 do j=jsc,jec ; do i=isc,iec
3477
3478 if ((iss%hmask(i,j) == 1) .OR. (iss%hmask(i,j) == 3)) then
3479
3480 if (trim(cs%ice_viscosity_compute) == "CONSTANT") then
3481 cs%ice_visc(i,j,1) = 1e15 * (us%kg_m3_to_R*us%m_to_L*us%m_s_to_L_T) * &
3482 (g%areaT(i,j) * max(iss%h_shelf(i,j),cs%min_h_shelf))
3483 ! constant viscocity for debugging
3484 elseif (trim(cs%ice_viscosity_compute) == "OBS") then
3485 if (cs%AGlen_visc(i,j) >0) then
3486 cs%ice_visc(i,j,1) = (g%areaT(i,j) * max(iss%h_shelf(i,j),cs%min_h_shelf)) * &
3487 max(cs%AGlen_visc(i,j) ,cs%min_ice_visc)
3488 endif
3489 ! Here CS%Aglen_visc(i,j) is the ice viscosity [R L2 T-1 ~> Pa s] computed from obs and read from a file
3490 elseif (model_qp1) then
3491 ! calculate viscosity at 1 cell-centered quadrature point per cell
3492
3493 visc_coef = (cs%AGlen_visc(i,j))**(-1./n_g)
3494 ! Units of Aglen_visc [Pa-(n_g) s-1]
3495
3496 ux = ((u_shlf(i-1,j-1) * cs%PhiC(1,i,j)) + &
3497 (u_shlf(i,j) * cs%PhiC(7,i,j))) + &
3498 ((u_shlf(i-1,j) * cs%PhiC(5,i,j)) + &
3499 (u_shlf(i,j-1) * cs%PhiC(3,i,j)))
3500
3501 vx = ((v_shlf(i-1,j-1) * cs%PhiC(1,i,j)) + &
3502 (v_shlf(i,j) * cs%PhiC(7,i,j))) + &
3503 ((v_shlf(i-1,j) * cs%PhiC(5,i,j)) + &
3504 (v_shlf(i,j-1) * cs%PhiC(3,i,j)))
3505
3506 uy = ((u_shlf(i-1,j-1) * cs%PhiC(2,i,j)) + &
3507 (u_shlf(i,j) * cs%PhiC(8,i,j))) + &
3508 ((u_shlf(i-1,j) * cs%PhiC(6,i,j)) + &
3509 (u_shlf(i,j-1) * cs%PhiC(4,i,j)))
3510
3511 vy = ((v_shlf(i-1,j-1) * cs%PhiC(2,i,j)) + &
3512 (v_shlf(i,j) * cs%PhiC(8,i,j))) + &
3513 ((v_shlf(i-1,j) * cs%PhiC(6,i,j)) + &
3514 (v_shlf(i,j-1) * cs%PhiC(4,i,j)))
3515
3516 cs%ice_visc(i,j,1) = (g%areaT(i,j) * max(iss%h_shelf(i,j),cs%min_h_shelf)) * &
3517 max(0.5 * visc_coef * &
3518 (us%s_to_T**2 * (((ux**2) + (vy**2)) + ((ux*vy) + 0.25*((uy+vx)**2)) + eps_min**2))**((1.-n_g)/(2.*n_g)) * &
3519 (us%Pa_to_RL2_T2*us%s_to_T),cs%min_ice_visc) ! Rescale after the fractional power law.
3520 ! Store Newton tangent stiffness data: strain rates and coefficient for Newton iterations.
3521 ! The Newton correction coefficient is (1/n-1)/2 * ice_visc / eps_e2,
3522 ! where eps_e2 = ux^2 + vy^2 + ux*vy + (uy+vx)^2/4 + eps_min^2 [T-2].
3523 ! It is zero where ice_visc is limited by min_ice_visc (viscosity is not smooth there).
3524 cs%newton_str_ux(i,j,1) = ux ; cs%newton_str_vy(i,j,1) = vy
3525 cs%newton_str_sh(i,j,1) = uy + vx
3526 cs%newton_visc_factor(i,j,1) = 0.0
3527 if (cs%ice_visc(i,j,1) > cs%min_ice_visc * (g%areaT(i,j) * max(iss%h_shelf(i,j),cs%min_h_shelf))) then
3528 cs%newton_visc_factor(i,j,1) = (0.5*(1./n_g - 1.) / &
3529 (((ux**2) + (vy**2)) + ((ux*vy) + 0.25*((uy+vx)**2)) + eps_min**2)) * &
3530 cs%ice_visc(i,j,1)
3531 endif
3532 elseif (model_qp4) then
3533 !calculate viscosity at 4 quadrature points per cell
3534
3535 visc_coef = (cs%AGlen_visc(i,j))**(-1./n_g)
3536
3537 do iq=1,2 ; do jq=1,2
3538
3539 ux = ((u_shlf(i-1,j-1) * cs%Phi(1,2*(jq-1)+iq,i,j)) + &
3540 (u_shlf(i,j) * cs%Phi(7,2*(jq-1)+iq,i,j))) + &
3541 ((u_shlf(i,j-1) * cs%Phi(3,2*(jq-1)+iq,i,j)) + &
3542 (u_shlf(i-1,j) * cs%Phi(5,2*(jq-1)+iq,i,j)))
3543
3544 vx = ((v_shlf(i-1,j-1) * cs%Phi(1,2*(jq-1)+iq,i,j)) + &
3545 (v_shlf(i,j) * cs%Phi(7,2*(jq-1)+iq,i,j))) + &
3546 ((v_shlf(i,j-1) * cs%Phi(3,2*(jq-1)+iq,i,j)) + &
3547 (v_shlf(i-1,j) * cs%Phi(5,2*(jq-1)+iq,i,j)))
3548
3549 uy = ((u_shlf(i-1,j-1) * cs%Phi(2,2*(jq-1)+iq,i,j)) + &
3550 (u_shlf(i,j) * cs%Phi(8,2*(jq-1)+iq,i,j))) + &
3551 ((u_shlf(i,j-1) * cs%Phi(4,2*(jq-1)+iq,i,j)) + &
3552 (u_shlf(i-1,j) * cs%Phi(6,2*(jq-1)+iq,i,j)))
3553
3554 vy = ((v_shlf(i-1,j-1) * cs%Phi(2,2*(jq-1)+iq,i,j)) + &
3555 (v_shlf(i,j) * cs%Phi(8,2*(jq-1)+iq,i,j))) + &
3556 ((v_shlf(i,j-1) * cs%Phi(4,2*(jq-1)+iq,i,j)) + &
3557 (v_shlf(i-1,j) * cs%Phi(6,2*(jq-1)+iq,i,j)))
3558
3559 cs%ice_visc(i,j,2*(jq-1)+iq) = (g%areaT(i,j) * max(iss%h_shelf(i,j),cs%min_h_shelf)) * &
3560 max(0.5 * visc_coef * &
3561 (us%s_to_T**2*(((ux**2) + (vy**2)) + ((ux*vy) + 0.25*((uy+vx)**2)) + eps_min**2))**((1.-n_g)/(2.*n_g)) * &
3562 (us%Pa_to_RL2_T2*us%s_to_T),cs%min_ice_visc) ! Rescale after the fractional power law.
3563 ! Store Newton tangent stiffness data at each quadrature point.
3564 cs%newton_str_ux(i,j,2*(jq-1)+iq) = ux ; cs%newton_str_vy(i,j,2*(jq-1)+iq) = vy
3565 cs%newton_str_sh(i,j,2*(jq-1)+iq) = uy + vx
3566 cs%newton_visc_factor(i,j,2*(jq-1)+iq) = 0.0
3567 if (cs%ice_visc(i,j,2*(jq-1)+iq) > &
3568 cs%min_ice_visc * (g%areaT(i,j) * max(iss%h_shelf(i,j),cs%min_h_shelf))) then
3569 cs%newton_visc_factor(i,j,2*(jq-1)+iq) = (0.5*(1./n_g - 1.) / &
3570 (((ux**2) + (vy**2)) + ((ux*vy) + 0.25*((uy+vx)**2)) + eps_min**2)) * &
3571 cs%ice_visc(i,j,2*(jq-1)+iq)
3572 endif
3573 enddo ; enddo
3574 endif
3575 endif
3576 enddo ; enddo
3577
3578end subroutine calc_shelf_visc
3579
3580
3581!> Update basal shear
3582subroutine calc_shelf_taub(CS, ISS, G, US, u_shlf, v_shlf)
3583 type(ice_shelf_dyn_cs), intent(inout) :: CS !< A pointer to the ice shelf control structure
3584 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
3585 !! the ice-shelf state
3586 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
3587 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
3588 real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
3589 intent(inout) :: u_shlf !< The zonal ice shelf velocity [L T-1 ~> m s-1].
3590 real, dimension(G%IsdB:G%IedB,G%JsdB:G%JedB), &
3591 intent(inout) :: v_shlf !< The meridional ice shelf velocity [L T-1 ~> m s-1].
3592
3593! also this subroutine updates the nonlinear part of the basal traction
3594
3595! this may be subject to change later... to make it "hybrid"
3596
3597 integer :: i, j, iscq, iecq, jscq, jecq, isd, jsd, ied, jed, iegq, jegq
3598 integer :: giec, gjec, gisc, gjsc, isc, jsc, iec, jec, is, js
3599 real :: umid, vmid ! Velocities [L T-1 ~> m s-1]
3600 real :: eps_min ! A minimal strain rate used in the Glens flow law expression [T-1 ~> s-1]
3601 real :: unorm ! The magnitude of the velocity in mks units for use with fractional powers [m s-1]
3602 real :: alpha ! Coulomb coefficient [nondim]
3603 real :: Hf !"floatation thickness" for Coulomb friction [Z ~> m]
3604 real :: fN ! Effective pressure (ice pressure - ocean pressure) for Coulomb friction [R Z L T-2 ~> Pa]
3605 real :: fB !for Coulomb Friction [(T L-1)^CS%CF_PostPeak ~> (s m-1)^CS%CF_PostPeak]
3606 real :: fBuq ! fB * unorm^CF_PostPeak, for Coulomb Newton correction [nondim]
3607 real :: unorm_code2 ! Squared velocity magnitude in code units [L2 T-2 ~> m2 s-2]
3608
3609 isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
3610 iscq = g%iscB ; iecq = g%iecB ; jscq = g%jscB ; jecq = g%jecB
3611 isd = g%isd ; jsd = g%jsd ; ied = g%ied ; jed = g%jed
3612 iegq = g%iegB ; jegq = g%jegB
3613 gisc = g%domain%nihalo+1 ; gjsc = g%domain%njhalo+1
3614 giec = g%domain%niglobal+gisc ; gjec = g%domain%njglobal+gjsc
3615 is = iscq - 1 ; js = jscq - 1
3616
3617 eps_min = cs%eps_glen_min
3618
3619 if (cs%CoulombFriction) then
3620 if (cs%CF_PostPeak/=1.0) THEN
3621 alpha = (cs%CF_PostPeak-1.0)**(cs%CF_PostPeak-1.0) / cs%CF_PostPeak**cs%CF_PostPeak ![nondim]
3622 else
3623 alpha = 1.0
3624 endif
3625 endif
3626
3627 do j=jsd+1,jed
3628 do i=isd+1,ied
3629 cs%newton_drag_coef(i,j) = 0.0
3630 if ((iss%hmask(i,j) == 1) .OR. (iss%hmask(i,j) == 3)) then
3631 umid = ((u_shlf(i,j) + u_shlf(i-1,j-1)) + (u_shlf(i,j-1) + u_shlf(i-1,j))) * 0.25
3632 vmid = ((v_shlf(i,j) + v_shlf(i-1,j-1)) + (v_shlf(i,j-1) + v_shlf(i-1,j))) * 0.25
3633 unorm_code2 = ((umid**2) + (vmid**2)) + (eps_min**2 * ((g%dxT(i,j)**2) + (g%dyT(i,j)**2)))
3634 unorm = us%L_T_to_m_s * sqrt( unorm_code2 )
3635
3636 !Coulomb friction (Schoof 2005, Gagliardini et al 2007)
3637 if (cs%CoulombFriction) then
3638 !Effective pressure
3639 hf = max((cs%density_ocean_avg/cs%density_ice) * cs%bed_elev(i,j), 0.0)
3640 fn = max((us%L_to_Z*(cs%density_ice * cs%g_Earth) * (max(iss%h_shelf(i,j),cs%min_h_shelf) - hf)), cs%CF_MinN)
3641 fb = alpha * (cs%C_basal_friction(i,j) / (cs%CF_Max * fn))**(cs%CF_PostPeak/cs%n_basal_fric)
3642 fbuq = fb * unorm**cs%CF_PostPeak
3643
3644 cs%basal_traction(i,j) = ((g%areaT(i,j) * cs%C_basal_friction(i,j)) * &
3645 (unorm**(cs%n_basal_fric-1.0) / (1.0 + fbuq)**(cs%n_basal_fric))) * &
3646 us%L_T_to_m_s ! Restore the scaling after the fractional power law.
3647 else
3648 !linear (CS%n_basal_fric=1) or "Weertman"/power-law (CS%n_basal_fric /= 1)
3649 fbuq = 0.0
3650 cs%basal_traction(i,j) = ((g%areaT(i,j) * cs%C_basal_friction(i,j)) * (unorm**(cs%n_basal_fric-1))) * &
3651 us%L_T_to_m_s ! Rescale after the fractional power law.
3652 endif
3653
3654 cs%basal_traction(i,j)=max(cs%basal_traction(i,j), cs%min_basal_traction * g%areaT(i,j))
3655
3656 ! Store Newton basal drag data for Newton tangent stiffness correction.
3657 ! newton_drag_coef = 2 * d(basal_trac)/d(|u|^2),
3658 ! where d(tau_b_i)/d(u_j) = basal_trac*delta_ij + newton_drag_coef*u_i*u_j
3659 ! This is the coefficient of the rank-1 correction u_i*(u.delta_u) to the Picard basal stiffness.
3660 ! For Weertman: newton_drag_coef = (m-1) * basal_trac/|u|^2
3661 ! For Coulomb: newton_drag_coef = basal_trac/|u|^2 * [(m-1) - m*q*fB*|u|^q/(1+fB*|u|^q)]
3662 cs%newton_umid(i,j) = umid
3663 cs%newton_vmid(i,j) = vmid
3664 ! unorm_code2: squared velocity magnitude in code units [L2 T-2], including regularization
3665 ! (same expression as inside the sqrt in unorm, without US%L_T_to_m_s factor)
3666 if (cs%CoulombFriction) then
3667 cs%newton_drag_coef(i,j) = (1.0 / max(unorm_code2, epsilon(unorm_code2))) * &
3668 cs%basal_traction(i,j) * ((cs%n_basal_fric - 1.) - &
3669 cs%n_basal_fric * cs%CF_PostPeak * fbuq / (1. + fbuq))
3670 else
3671 cs%newton_drag_coef(i,j) = real(cs%n_basal_fric - 1.) * cs%basal_traction(i,j) / &
3672 max(unorm_code2, epsilon(unorm_code2))
3673 endif
3674 endif
3675 enddo
3676 enddo
3677
3678end subroutine calc_shelf_taub
3679
3680subroutine update_od_ffrac(CS, G, US, ocean_mass, find_avg)
3681 type(ice_shelf_dyn_cs), intent(inout) :: CS !< A pointer to the ice shelf control structure
3682 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
3683 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
3684 real, dimension(SZDI_(G),SZDJ_(G)), &
3685 intent(in) :: ocean_mass !< The mass per unit area of the ocean [R Z ~> kg m-2].
3686 logical, intent(in) :: find_avg !< If true, find the average of OD and ffrac, and
3687 !! reset the underlying running sums to 0.
3688
3689 integer :: isc, iec, jsc, jec, i, j
3690 real :: I_rho_ocean ! A typical specific volume of the ocean [R-1 ~> m3 kg-1]
3691 real :: I_counter
3692
3693 i_rho_ocean = 1.0 / cs%density_ocean_avg
3694
3695 isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
3696
3697 do j=jsc,jec ; do i=isc,iec
3698 cs%OD_rt(i,j) = cs%OD_rt(i,j) + ocean_mass(i,j)*i_rho_ocean
3699 if (ocean_mass(i,j)*i_rho_ocean > cs%thresh_float_col_depth) then
3700 cs%ground_frac_rt(i,j) = cs%ground_frac_rt(i,j) + 1.0
3701 endif
3702 enddo ; enddo
3703 cs%OD_rt_counter = cs%OD_rt_counter + 1
3704
3705 if (find_avg) then
3706 i_counter = 1.0 / real(cs%OD_rt_counter)
3707 do j=jsc,jec ; do i=isc,iec
3708 cs%ground_frac(i,j) = 1.0 - (cs%ground_frac_rt(i,j) * i_counter)
3709 cs%OD_av(i,j) = cs%OD_rt(i,j) * i_counter
3710
3711 cs%OD_rt(i,j) = 0.0 ; cs%ground_frac_rt(i,j) = 0.0 ; cs%OD_rt_counter = 0
3712 enddo ; enddo
3713
3714 call pass_var(cs%ground_frac, g%domain, complete=.false.)
3715 call pass_var(cs%OD_av, g%domain, complete=.true.)
3716 endif
3717
3718end subroutine update_od_ffrac
3719
3720subroutine update_od_ffrac_uncoupled(CS, G, h_shelf)
3721 type(ice_shelf_dyn_cs), intent(inout) :: CS !< A pointer to the ice shelf control structure
3722 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
3723 real, dimension(SZDI_(G),SZDJ_(G)), &
3724 intent(in) :: h_shelf !< the thickness of the ice shelf [Z ~> m].
3725
3726 integer :: i, j, isd, ied, jsd, jed
3727 real :: rhoi_rhow, OD
3728
3729 rhoi_rhow = cs%density_ice / cs%density_ocean_avg
3730 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
3731
3732 do j=jsd,jed
3733 do i=isd,ied
3734 od = cs%bed_elev(i,j) - rhoi_rhow * max(h_shelf(i,j),cs%min_h_shelf)
3735 if (od >= 0) then
3736 ! ice thickness does not take up whole ocean column -> floating
3737 cs%OD_av(i,j) = od
3738 cs%ground_frac(i,j) = 0.
3739 else
3740 cs%OD_av(i,j) = 0.
3741 cs%ground_frac(i,j) = 1.
3742 endif
3743 enddo
3744 enddo
3745
3746end subroutine update_od_ffrac_uncoupled
3747
3748subroutine change_in_draft(CS, G, h_shelf0, h_shelf1, ddraft)
3749 type(ice_shelf_dyn_cs), intent(inout) :: cs !< A pointer to the ice shelf control structure
3750 type(ocean_grid_type), intent(in) :: g !< The grid structure used by the ice shelf.
3751 real, dimension(SZDI_(G),SZDJ_(G)), &
3752 intent(in) :: h_shelf0 !< the previous thickness of the ice shelf [Z ~> m].
3753 real, dimension(SZDI_(G),SZDJ_(G)), &
3754 intent(in) :: h_shelf1 !< the current thickness of the ice shelf [Z ~> m].
3755 real, dimension(SZDI_(G),SZDJ_(G)), &
3756 intent(inout) :: ddraft !< the change in shelf draft thickness
3757 real :: b0,b1
3758 integer :: i, j, isc, iec, jsc, jec
3759 real :: rhoi_rhow, od
3760
3761 rhoi_rhow = cs%density_ice / cs%density_ocean_avg
3762 isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
3763 ddraft = 0.0
3764
3765 do j=jsc,jec
3766 do i=isc,iec
3767
3768 b0 = 0.0 ; b1 = 0.0
3769
3770 if (h_shelf0(i,j)>0.0) then
3771 od = cs%bed_elev(i,j) - rhoi_rhow * h_shelf0(i,j)
3772 if (od >= 0) then
3773 !floating
3774 b0 = rhoi_rhow * h_shelf0(i,j)
3775 else
3776 b0 = cs%bed_elev(i,j)
3777 endif
3778 endif
3779
3780 if (h_shelf1(i,j)>0.0) then
3781 od = cs%bed_elev(i,j) - rhoi_rhow * h_shelf1(i,j)
3782 if (od >= 0) then
3783 !floating
3784 b1 = rhoi_rhow * h_shelf1(i,j)
3785 else
3786 b1 = cs%bed_elev(i,j)
3787 endif
3788 endif
3789
3790 ddraft(i,j) = b1-b0
3791 enddo
3792 enddo
3793end subroutine change_in_draft
3794
3795!> This subroutine calculates the gradients of bilinear basis elements that
3796!! that are centered at the vertices of the cell. Values are calculated at
3797!! points of gaussian quadrature.
3798subroutine bilinear_shape_functions (X, Y, Phi, area)
3799 real, dimension(4), intent(in) :: X !< The x-positions of the vertices of the quadrilateral [L ~> m].
3800 real, dimension(4), intent(in) :: Y !< The y-positions of the vertices of the quadrilateral [L ~> m].
3801 real, dimension(8,4), intent(inout) :: Phi !< The gradients of bilinear basis elements at Gaussian
3802 !! quadrature points surrounding the cell vertices [L-1 ~> m-1].
3803 real, intent(out) :: area !< The quadrilateral cell area [L2 ~> m2].
3804
3805! X and Y must be passed in the form
3806 ! 3 - 4
3807 ! | |
3808 ! 1 - 2
3809
3810! this subroutine calculates the gradients of bilinear basis elements that
3811! that are centered at the vertices of the cell. values are calculated at
3812! points of gaussian quadrature. (in 1D: .5 * (1 +/- sqrt(1/3)) for [0,1])
3813! (ordered in same way as vertices)
3814!
3815! Phi(2*i-1,j) gives d(Phi_i)/dx at quadrature point j
3816! Phi(2*i,j) gives d(Phi_i)/dy at quadrature point j
3817! Phi_i is equal to 1 at vertex i, and 0 at vertex k /= i, and bilinear
3818!
3819! This should be a one-off; once per nonlinear solve? once per lifetime?
3820! ... will all cells have the same shape and dimension?
3821
3822 real, dimension(4) :: xquad, yquad ! [nondim]
3823 real :: a,b,c,d ! Various lengths [L ~> m]
3824 real :: xexp, yexp ! [nondim]
3825 integer :: node, qpoint, xnode, ynode
3826
3827 xquad(1:3:2) = .5 * (1-sqrt(1./3)) ; yquad(1:2) = .5 * (1-sqrt(1./3))
3828 xquad(2:4:2) = .5 * (1+sqrt(1./3)) ; yquad(3:4) = .5 * (1+sqrt(1./3))
3829
3830 do qpoint=1,4
3831
3832 a = ((-x(1)*(1-yquad(qpoint)))+(x(4)*yquad(qpoint))) + ((x(2)*(1-yquad(qpoint)))-(x(3)*yquad(qpoint))) !d(x)/d(x*)
3833 b = ((-y(1)*(1-yquad(qpoint)))+(y(4)*yquad(qpoint))) + ((y(2)*(1-yquad(qpoint)))-(y(3)*yquad(qpoint))) !d(y)/d(x*)
3834 c = ((-x(1)*(1-xquad(qpoint)))+(x(4)*xquad(qpoint))) + ((-x(2)*xquad(qpoint))+(x(3)*(1-xquad(qpoint))))!d(x)/d(y*)
3835 d = ((-y(1)*(1-xquad(qpoint)))+(y(4)*xquad(qpoint))) + ((-y(2)*xquad(qpoint))+(y(3)*(1-xquad(qpoint))))!d(y)/d(y*)
3836
3837 do node=1,4
3838
3839 xnode = 2-mod(node,2) ; ynode = ceiling(real(node)/2)
3840
3841 if (ynode == 1) then
3842 yexp = 1-yquad(qpoint)
3843 else
3844 yexp = yquad(qpoint)
3845 endif
3846
3847 if (1 == xnode) then
3848 xexp = 1-xquad(qpoint)
3849 else
3850 xexp = xquad(qpoint)
3851 endif
3852
3853 phi(2*node-1,qpoint) = ( d * (2 * xnode - 3) * yexp - b * (2 * ynode - 3) * xexp) / ((a*d)-(b*c))
3854 phi(2*node,qpoint) = (-c * (2 * xnode - 3) * yexp + a * (2 * ynode - 3) * xexp) / ((a*d)-(b*c))
3855
3856 enddo
3857 enddo
3858
3859 area = quad_area(x, y)
3860
3861end subroutine bilinear_shape_functions
3862
3863!> This subroutine calculates the gradients of bilinear basis elements that are centered at the
3864!! vertices of the cell using a locally orthogoal MOM6 grid. Values are calculated at
3865!! points of gaussian quadrature.
3866subroutine bilinear_shape_fn_grid(G, i, j, Phi)
3867 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
3868 integer, intent(in) :: i !< The i-index in the grid to work on.
3869 integer, intent(in) :: j !< The j-index in the grid to work on.
3870 real, dimension(8,4), intent(inout) :: Phi !< The gradients of bilinear basis elements at Gaussian
3871 !! quadrature points surrounding the cell vertices [L-1 ~> m-1].
3872
3873! This subroutine calculates the gradients of bilinear basis elements that
3874! that are centered at the vertices of the cell. The values are calculated at
3875! points of gaussian quadrature. (in 1D: .5 * (1 +/- sqrt(1/3)) for [0,1])
3876! (ordered in same way as vertices)
3877!
3878! Phi(2*i-1,j) gives d(Phi_i)/dx at quadrature point j
3879! Phi(2*i,j) gives d(Phi_i)/dy at quadrature point j
3880! Phi_i is equal to 1 at vertex i, and 0 at vertex k /= i, and bilinear
3881!
3882! This should be a one-off; once per nonlinear solve? once per lifetime?
3883
3884 real, dimension(4) :: xquad, yquad ! [nondim]
3885 real :: a, d ! Interpolated grid spacings [L ~> m]
3886 real :: xexp, yexp ! [nondim]
3887 integer :: node, qpoint, xnode, ynode
3888
3889 xquad(1:3:2) = .5 * (1-sqrt(1./3)) ; yquad(1:2) = .5 * (1-sqrt(1./3))
3890 xquad(2:4:2) = .5 * (1+sqrt(1./3)) ; yquad(3:4) = .5 * (1+sqrt(1./3))
3891
3892 do qpoint=1,4
3893 if (j>1) then
3894 a = (g%dxCv(i,j-1) * (1-yquad(qpoint))) + (g%dxCv(i,j) * yquad(qpoint)) ! d(x)/d(x*)
3895 else
3896 a = g%dxCv(i,j) !* yquad(qpoint) ! d(x)/d(x*)
3897 endif
3898 if (i>1) then
3899 d = (g%dyCu(i-1,j) * (1-xquad(qpoint))) + (g%dyCu(i,j) * xquad(qpoint)) ! d(y)/d(y*)
3900 else
3901 d = g%dyCu(i,j) !* xquad(qpoint)
3902 endif
3903! a = G%dxCv(i,J-1) * (1-yquad(qpoint)) + G%dxCv(i,J) * yquad(qpoint) ! d(x)/d(x*)
3904! d = G%dyCu(I-1,j) * (1-xquad(qpoint)) + G%dyCu(I,j) * xquad(qpoint) ! d(y)/d(y*)
3905
3906 do node=1,4
3907 xnode = 2-mod(node,2) ; ynode = ceiling(real(node)/2)
3908
3909 if (ynode == 1) then
3910 yexp = 1-yquad(qpoint)
3911 else
3912 yexp = yquad(qpoint)
3913 endif
3914
3915 if (1 == xnode) then
3916 xexp = 1-xquad(qpoint)
3917 else
3918 xexp = xquad(qpoint)
3919 endif
3920
3921 phi(2*node-1,qpoint) = ( (d * (2 * xnode - 3)) * yexp ) / (a*d)
3922 phi(2*node,qpoint) = ( (a * (2 * ynode - 3)) * xexp ) / (a*d)
3923
3924 enddo
3925 enddo
3926
3927end subroutine bilinear_shape_fn_grid
3928
3929!> This subroutine calculates the gradients of bilinear basis elements that are centered at the
3930!! vertices of the cell using a locally orthogoal MOM6 grid. Values are calculated at
3931!! a sinlge cell-centered quadrature point, which should match the grid cell h-point
3932subroutine bilinear_shape_fn_grid_1qp(G, i, j, Phi)
3933 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
3934 integer, intent(in) :: i !< The i-index in the grid to work on.
3935 integer, intent(in) :: j !< The j-index in the grid to work on.
3936 real, dimension(8), intent(inout) :: Phi !< The gradients of bilinear basis elements at Gaussian
3937 !! quadrature points surrounding the cell vertices [L-1 ~> m-1].
3938
3939! This subroutine calculates the gradients of bilinear basis elements that
3940! that are centered at the vertices of the cell. The values are calculated at
3941! a cell-cented point of gaussian quadrature. (in 1D: .5 for [0,1])
3942! (ordered in same way as vertices)
3943!
3944! Phi(2*i-1) gives d(Phi_i)/dx at the quadrature point
3945! Phi(2*i) gives d(Phi_i)/dy at the quadrature point
3946! Phi_i is equal to 1 at vertex i, and 0 at vertex k /= i, and bilinear
3947
3948 real :: a, d ! Interpolated grid spacings [L ~> m]
3949 real :: xexp=0.5, yexp=0.5 ! [nondim]
3950 integer :: node, qpoint, xnode, ynode
3951
3952 ! d(x)/d(x*)
3953 if (j>1) then
3954 a = 0.5 * (g%dxCv(i,j-1) + g%dxCv(i,j))
3955 else
3956 a = g%dxCv(i,j)
3957 endif
3958
3959 ! d(y)/d(y*)
3960 if (i>1) then
3961 d = 0.5 * (g%dyCu(i-1,j) + g%dyCu(i,j))
3962 else
3963 d = g%dyCu(i,j)
3964 endif
3965
3966 do node=1,4
3967 xnode = 2-mod(node,2) ; ynode = ceiling(real(node)/2)
3968 phi(2*node-1) = ( (d * (2 * xnode - 3)) * yexp ) / (a*d)
3969 phi(2*node) = ( (a * (2 * ynode - 3)) * xexp ) / (a*d)
3970 enddo
3971end subroutine bilinear_shape_fn_grid_1qp
3972
3973
3974subroutine bilinear_shape_functions_subgrid(Phisub, nsub)
3975 integer, intent(in) :: nsub !< The number of subgridscale quadrature locations in each direction
3976 real, dimension(2,2,nsub,nsub,2,2), &
3977 intent(inout) :: Phisub !< Quadrature structure weights at subgridscale
3978 !! locations for finite element calculations [nondim]
3979
3980 ! this subroutine is a helper for interpolation of floatation condition
3981 ! for the purposes of evaluating the terms \int (u,v) \phi_i dx dy in a cell that is
3982 ! in partial floatation
3983 ! the array Phisub contains the values of \phi_i (where i is a node of the cell)
3984 ! at quad point j
3985 ! i think this general approach may not work for nonrectangular elements...
3986 !
3987
3988 ! Phisub(q1,q2,i,j,k,l)
3989 ! q1: quad point x-index
3990 ! q2: quad point y-index
3991 ! i: subgrid index in x-direction
3992 ! j: subgrid index in y-direction
3993 ! k: basis function x-index
3994 ! l: basis function y-index
3995
3996 ! e.g. k=1,l=1 => node 1
3997 ! q1=2,q2=1 => quad point 2
3998
3999 ! 3 - 4
4000 ! | |
4001 ! 1 - 2
4002
4003 integer :: i, j, qx, qy
4004 real,dimension(2) :: xquad
4005 real :: x0, y0, x, y, fracx
4006
4007 xquad(1) = .5 * (1-sqrt(1./3)) ; xquad(2) = .5 * (1+sqrt(1./3))
4008 fracx = 1.0/real(nsub)
4009
4010 do j=1,nsub ; do i=1,nsub
4011 x0 = (i-1) * fracx ; y0 = (j-1) * fracx
4012 do qy=1,2 ; do qx=1,2
4013 x = x0 + fracx*xquad(qx)
4014 y = y0 + fracx*xquad(qy)
4015 phisub(qx,qy,i,j,1,1) = (1.0-x) * (1.0-y)
4016 phisub(qx,qy,i,j,1,2) = (1.0-x) * y
4017 phisub(qx,qy,i,j,2,1) = x * (1.0-y)
4018 phisub(qx,qy,i,j,2,2) = x * y
4019 enddo ; enddo
4020 enddo ; enddo
4021
4022end subroutine bilinear_shape_functions_subgrid
4023
4024
4025subroutine update_velocity_masks(CS, G, hmask, umask, vmask, u_face_mask, v_face_mask)
4026 type(ice_shelf_dyn_cs),intent(in) :: CS !< A pointer to the ice shelf dynamics control structure
4027 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
4028 real, dimension(SZDI_(G),SZDJ_(G)), &
4029 intent(in) :: hmask !< A mask indicating which tracer points are
4030 !! partly or fully covered by an ice-shelf
4031 real, dimension(SZDIB_(G),SZDJB_(G)), &
4032 intent(out) :: umask !< A coded mask indicating the nature of the
4033 !! zonal flow at the corner point
4034 real, dimension(SZDIB_(G),SZDJB_(G)), &
4035 intent(out) :: vmask !< A coded mask indicating the nature of the
4036 !! meridional flow at the corner point
4037 real, dimension(SZDIB_(G),SZDJB_(G)), &
4038 intent(out) :: u_face_mask !< A coded mask for velocities at the C-grid u-face
4039 real, dimension(SZDIB_(G),SZDJB_(G)), &
4040 intent(out) :: v_face_mask !< A coded mask for velocities at the C-grid v-face
4041 ! sets masks for velocity solve
4042 ! ignores the fact that their might be ice-free cells - this only considers the computational boundary
4043
4044 ! !!!IMPORTANT!!! relies on thickness mask - assumed that this is called after hmask has been updated & halo-updated
4045
4046 integer :: i, j, k, iscq, iecq, jscq, jecq, isd, jsd, is, js, iegq, jegq
4047 integer :: giec, gjec, gisc, gjsc, isc, jsc, iec, jec
4048
4049 isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
4050 iscq = g%iscB ; iecq = g%iecB ; jscq = g%jscB ; jecq = g%jecB
4051 isd = g%isd ; jsd = g%jsd
4052 iegq = g%iegB ; jegq = g%jegB
4053 gisc = g%Domain%nihalo ; gjsc = g%Domain%njhalo
4054 giec = g%Domain%niglobal+gisc ; gjec = g%Domain%njglobal+gjsc
4055
4056 umask(:,:) = 0 ; vmask(:,:) = 0
4057 u_face_mask(:,:) = 0 ; v_face_mask(:,:) = 0
4058
4059 if (g%symmetric) then
4060 is = isd ; js = jsd
4061 else
4062 is = isd+1 ; js = jsd+1
4063 endif
4064
4065 do j=js,g%jed ; do i=is,g%ied
4066 if (hmask(i,j) == 1 .or. hmask(i,j)==3) then
4067 umask(i-1:i,j-1:j)=1
4068 vmask(i-1:i,j-1:j)=1
4069 endif
4070 enddo ; enddo
4071
4072 do j=js,g%jed
4073 do i=is,g%ied
4074
4075 if ((hmask(i,j) == 1) .OR. (hmask(i,j) == 3)) then
4076
4077 do k=0,1
4078
4079 select case (int(cs%u_face_mask_bdry(i-1+k,j)))
4080 case (5)
4081 umask(i-1+k,j-1:j) = 3.
4082 u_face_mask(i-1+k,j) = 5.
4083 case (3)
4084 umask(i-1+k,j-1:j) = 3.
4085 vmask(i-1+k,j-1:j) = 3.
4086 u_face_mask(i-1+k,j) = 3.
4087 case (6)
4088 vmask(i-1+k,j-1:j) = 3.
4089 u_face_mask(i-1+k,j) = 6.
4090 case (2)
4091 u_face_mask(i-1+k,j) = 2.
4092 case (4)
4093 umask(i-1+k,j-1:j) = 0.
4094 u_face_mask(i-1+k,j) = 4.
4095 case (0)
4096 umask(i-1+k,j-1:j) = 0.
4097 u_face_mask(i-1+k,j) = 0.
4098 case (1) ! stress free x-boundary
4099 umask(i-1+k,j-1:j) = 0.
4100 case default
4101 umask(i-1+k,j-1) = max(1. , umask(i-1+k,j-1))
4102 umask(i-1+k,j) = max(1. , umask(i-1+k,j))
4103 end select
4104 enddo
4105
4106 do k=0,1
4107
4108 select case (int(cs%v_face_mask_bdry(i,j-1+k)))
4109 case (5)
4110 vmask(i-1:i,j-1+k) = 3.
4111 v_face_mask(i,j-1+k) = 5.
4112 case (3)
4113 vmask(i-1:i,j-1+k) = 3.
4114 umask(i-1:i,j-1+k) = 3.
4115 v_face_mask(i,j-1+k) = 3.
4116 case (6)
4117 umask(i-1:i,j-1+k) = 3.
4118 v_face_mask(i,j-1+k) = 6.
4119 case (2)
4120 v_face_mask(i,j-1+k) = 2.
4121 case (4)
4122 vmask(i-1:i,j-1+k) = 0.
4123 v_face_mask(i,j-1+k) = 4.
4124 case (0)
4125 vmask(i-1:i,j-1+k) = 0.
4126 v_face_mask(i,j-1+k) = 0.
4127 case (1) ! stress free y-boundary
4128 vmask(i-1:i,j-1+k) = 0.
4129 case default
4130 vmask(i-1,j-1+k) = max(1. , vmask(i-1,j-1+k))
4131 vmask(i,j-1+k) = max(1. , vmask(i,j-1+k))
4132 end select
4133 enddo
4134
4135
4136 if (i < g%ied) then
4137 if ((hmask(i+1,j) == 0) .OR. (hmask(i+1,j) == 2)) then
4138 ! east boundary or adjacent to unfilled cell
4139 u_face_mask(i,j) = 2.
4140 endif
4141 endif
4142
4143 if (i > g%isd) then
4144 if ((hmask(i-1,j) == 0) .OR. (hmask(i-1,j) == 2)) then
4145 !adjacent to unfilled cell
4146 u_face_mask(i-1,j) = 2.
4147 endif
4148 endif
4149
4150 if (j > g%jsd) then
4151 if ((hmask(i,j-1) == 0) .OR. (hmask(i,j-1) == 2)) then
4152 !adjacent to unfilled cell
4153 v_face_mask(i,j-1) = 2.
4154 endif
4155 endif
4156
4157 if (j < g%jed) then
4158 if ((hmask(i,j+1) == 0) .OR. (hmask(i,j+1) == 2)) then
4159 !adjacent to unfilled cell
4160 v_face_mask(i,j) = 2.
4161 endif
4162 endif
4163
4164
4165 endif
4166
4167 enddo
4168 enddo
4169
4170 ! note: if the grid is nonsymmetric, there is a part that will not be transferred with a halo update
4171 ! so this subroutine must update its own symmetric part of the halo
4172
4173 call pass_vector(u_face_mask, v_face_mask, g%domain, to_all, cgrid_ne)
4174 call pass_vector(umask, vmask, g%domain, to_all, bgrid_ne)
4175
4176end subroutine update_velocity_masks
4177
4178!> Interpolate the ice shelf thickness from tracer point to nodal points,
4179!! subject to a mask.
4180subroutine interpolate_h_to_b(G, h_shelf, hmask, H_node, min_h_shelf)
4181 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
4182 real, dimension(SZDI_(G),SZDJ_(G)), &
4183 intent(in) :: h_shelf !< The ice shelf thickness at tracer points [Z ~> m].
4184 real, dimension(SZDI_(G),SZDJ_(G)), &
4185 intent(in) :: hmask !< A mask indicating which tracer points are
4186 !! partly or fully covered by an ice-shelf
4187 real, dimension(SZDIB_(G),SZDJB_(G)), &
4188 intent(inout) :: H_node !< The ice shelf thickness at nodal (corner)
4189 !! points [Z ~> m].
4190 real, intent(in) :: min_h_shelf !< The minimum ice thickness used during ice dynamics [Z ~> m].
4191
4192 integer :: i, j, isc, iec, jsc, jec, num_h, k, l, ic, jc
4193 real :: h_arr(2,2)
4194
4195 isc = g%isc ; jsc = g%jsc ; iec = g%iec ; jec = g%jec
4196
4197 h_node(:,:) = 0.0
4198
4199 ! H_node is node-centered; average over all cells that share that node
4200 ! if no (active) cells share the node then its value there is irrelevant
4201
4202 do j=jsc-1,jec
4203 do i=isc-1,iec
4204 num_h = 0
4205 do l=1,2 ; jc=j-1+l ; do k=1,2 ; ic=i-1+k
4206 if (hmask(ic,jc) == 1.0 .or. hmask(ic,jc) == 3.0) then
4207 h_arr(k,l)=max(h_shelf(ic,jc),min_h_shelf)
4208 num_h = num_h + 1
4209 else
4210 h_arr(k,l)=0.0
4211 endif
4212 if (num_h > 0) then
4213 h_node(i,j) = ((h_arr(1,1)+h_arr(2,2))+(h_arr(1,2)+h_arr(2,1))) / num_h
4214 endif
4215 enddo ; enddo
4216 enddo
4217 enddo
4218
4219 call pass_var(h_node, g%domain,position=corner)
4220
4221end subroutine interpolate_h_to_b
4222
4223!> Deallocates all memory associated with the ice shelf dynamics module
4224subroutine ice_shelf_dyn_end(CS)
4225 type(ice_shelf_dyn_cs), pointer :: cs !< A pointer to the ice shelf dynamics control structure
4226
4227 if (.not.associated(cs)) return
4228
4229 deallocate(cs%u_shelf, cs%v_shelf)
4230 deallocate(cs%taudx_shelf, cs%taudy_shelf)
4231 deallocate(cs%t_shelf, cs%tmask)
4232 deallocate(cs%u_bdry_val, cs%v_bdry_val)
4233 deallocate(cs%u_face_mask, cs%v_face_mask)
4234 deallocate(cs%umask, cs%vmask)
4235 deallocate(cs%u_face_mask_bdry, cs%v_face_mask_bdry)
4236 deallocate(cs%h_bdry_val)
4237 deallocate(cs%float_cond)
4238
4239 deallocate(cs%ice_visc, cs%AGlen_visc)
4240 deallocate(cs%newton_visc_factor, cs%newton_str_ux, cs%newton_str_vy, cs%newton_str_sh)
4241 deallocate(cs%newton_umid, cs%newton_vmid, cs%newton_drag_coef)
4242 deallocate(cs%basal_traction,cs%C_basal_friction)
4243 deallocate(cs%OD_rt, cs%OD_av)
4244 deallocate(cs%t_bdry_val, cs%bed_elev)
4245 deallocate(cs%ground_frac, cs%ground_frac_rt)
4246
4247 deallocate(cs)
4248
4249end subroutine ice_shelf_dyn_end
4250
4251
4252!> This subroutine updates the vertically averaged ice shelf temperature.
4253subroutine ice_shelf_temp(CS, ISS, G, US, time_step, melt_rate, Time)
4254 type(ice_shelf_dyn_cs), intent(inout) :: CS !< A pointer to the ice shelf control structure
4255 type(ice_shelf_state), intent(in) :: ISS !< A structure with elements that describe
4256 !! the ice-shelf state
4257 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
4258 type(unit_scale_type), intent(in) :: US !< A structure containing unit conversion factors
4259 real, intent(in) :: time_step !< The time step for this update [T ~> s].
4260 real, dimension(SZDI_(G),SZDJ_(G)), &
4261 intent(in) :: melt_rate !< basal melt rate [R Z T-1 ~> kg m-2 s-1]
4262 type(time_type), intent(in) :: Time !< The current model time
4263
4264! This subroutine takes the velocity (on the Bgrid) and timesteps
4265! (HT)_t = - div (uHT) + (adot Tsurf -bdot Tbot) once and then calculates T=HT/H
4266!
4267! The flux overflows are included here. That is because they will be used to advect 3D scalars
4268! into partial cells
4269
4270 real, dimension(SZDI_(G),SZDJ_(G)) :: th_after_uflux, th_after_vflux, TH ! Integrated temperatures [C Z ~> degC m]
4271 integer :: isd, ied, jsd, jed, i, j, isc, iec, jsc, jec
4272 real :: Tsurf ! Surface air temperature [C ~> degC]. This is hard coded but should be an input argument.
4273 real :: adot ! A surface heat exchange coefficient [R Z T-1 ~> kg m-2 s-1].
4274
4275
4276 ! For now adot and Tsurf are defined here adot=surf acc 0.1m/yr, Tsurf=-20oC, vary them later
4277 adot = (0.1/(365.0*86400.0))*us%m_to_Z*us%T_to_s * cs%density_ice
4278 tsurf = -20.0*us%degC_to_C
4279
4280 isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
4281 isc = g%isc ; iec = g%iec ; jsc = g%jsc ; jec = g%jec
4282
4283 th_after_uflux(:,:) = 0.0
4284 th_after_vflux(:,:) = 0.0
4285
4286 do j=jsd,jed ; do i=isd,ied
4287! if (ISS%hmask(i,j) > 1) then
4288 if ((iss%hmask(i,j) == 3) .or. (iss%hmask(i,j) == -2)) then
4289 cs%t_shelf(i,j) = cs%t_bdry_val(i,j)
4290 endif
4291 enddo ; enddo
4292
4293 do j=jsd,jed ; do i=isd,ied
4294 ! Convert the averge temperature to a depth integrated temperature.
4295 th(i,j) = cs%t_shelf(i,j)*iss%h_shelf(i,j)
4296 enddo ; enddo
4297
4298
4299 call ice_shelf_advect_temp_x(cs, g, time_step, iss%hmask, th, th_after_uflux)
4300 call ice_shelf_advect_temp_y(cs, g, time_step, iss%hmask, th_after_uflux, th_after_vflux)
4301
4302 do j=jsc,jec ; do i=isc,iec
4303 ! Convert the integrated temperature back to the average temperature.
4304! if ((ISS%hmask(i,j) == 1) .or. (ISS%hmask(i,j) == 2)) then
4305 if (iss%h_shelf(i,j) > 0.0) then
4306 cs%t_shelf(i,j) = th_after_vflux(i,j) / iss%h_shelf(i,j)
4307 else
4308 cs%t_shelf(i,j) = cs%T_shelf_missing
4309 endif
4310! endif
4311
4312 if ((iss%hmask(i,j) == 1) .or. (iss%hmask(i,j) == 2)) then
4313 if (iss%h_shelf(i,j) > 0.0) then
4314 cs%t_shelf(i,j) = cs%t_shelf(i,j) + &
4315 time_step*(adot*tsurf - melt_rate(i,j)*iss%tfreeze(i,j))/(cs%density_ice*iss%h_shelf(i,j))
4316 else
4317 ! the ice is about to melt away in this case set thickness, area, and mask to zero
4318 ! NOTE: not mass conservative, should maybe scale salt & heat flux for this cell
4319 cs%t_shelf(i,j) = cs%T_shelf_missing
4320 cs%tmask(i,j) = 0.0
4321 endif
4322 elseif (iss%hmask(i,j) == 0) then
4323 cs%t_shelf(i,j) = cs%T_shelf_missing
4324 elseif ((iss%hmask(i,j) == 3) .or. (iss%hmask(i,j) == -2)) then
4325 cs%t_shelf(i,j) = cs%t_bdry_val(i,j)
4326 endif
4327 enddo ; enddo
4328
4329 call pass_var(cs%t_shelf, g%domain, complete=.false.)
4330 call pass_var(cs%tmask, g%domain, complete=.true.)
4331
4332 if (cs%debug) then
4333 call hchksum(cs%t_shelf, "temp after front", g%HI, haloshift=3, unscale=us%C_to_degC)
4334 endif
4335
4336end subroutine ice_shelf_temp
4337
4338
4339subroutine ice_shelf_advect_temp_x(CS, G, time_step, hmask, h0, h_after_uflux)
4340 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
4341 type(ocean_grid_type), intent(inout) :: G !< The grid structure used by the ice shelf.
4342 real, intent(in) :: time_step !< The time step for this update [T ~> s].
4343 real, dimension(SZDI_(G),SZDJ_(G)), &
4344 intent(in) :: hmask !< A mask indicating which tracer points are
4345 !! partly or fully covered by an ice-shelf
4346 real, dimension(SZDI_(G),SZDJ_(G)), &
4347 intent(in) :: h0 !< The initial ice shelf thicknesses times temperature [C Z ~> degC m]
4348 real, dimension(SZDI_(G),SZDJ_(G)), &
4349 intent(inout) :: h_after_uflux !< The ice shelf thicknesses times temperature after
4350 !! the zonal mass fluxes [C Z ~> degC m]
4351
4352 ! use will be made of ISS%hmask here - its value at the boundary will be zero, just like uncovered cells
4353 ! if there is an input bdry condition, the thickness there will be set in initialization
4354
4355 integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
4356 integer :: i_off, j_off
4357 logical :: at_east_bdry, at_west_bdry
4358 real, dimension(-2:2) :: stencil ! A copy of the neighboring thicknesses times temperatures [C Z ~> degC m]
4359 real :: u_face ! Zonal velocity at a face, positive if out [L T-1 ~> m s-1]
4360 real :: flux_diff ! The difference in fluxes [C Z ~> degC m]
4361 real :: phi ! A limiting ratio [nondim]
4362
4363 is = g%isc-2 ; ie = g%iec+2 ; js = g%jsc ; je = g%jec ; isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
4364 i_off = g%idg_offset ; j_off = g%jdg_offset
4365
4366 do j=jsd+1,jed-1
4367 if (((j+j_off) <= g%domain%njglobal+g%domain%njhalo) .AND. &
4368 ((j+j_off) >= g%domain%njhalo+1)) then ! based on mehmet's code - only if btw north & south boundaries
4369
4370 stencil(:) = 0.0 ! This is probably unnecessary, as the code is written
4371! if (i+i_off == G%domain%nihalo+G%domain%nihalo)
4372 do i=is,ie
4373
4374 if (((i+i_off) <= g%domain%niglobal+g%domain%nihalo) .AND. &
4375 ((i+i_off) >= g%domain%nihalo+1)) then
4376
4377 if (i+i_off == g%domain%nihalo+1) then
4378 at_west_bdry=.true.
4379 else
4380 at_west_bdry=.false.
4381 endif
4382
4383 if (i+i_off == g%domain%niglobal+g%domain%nihalo) then
4384 at_east_bdry=.true.
4385 else
4386 at_east_bdry=.false.
4387 endif
4388
4389 if (hmask(i,j) == 1) then
4390
4391 h_after_uflux(i,j) = h0(i,j)
4392
4393 stencil(:) = h0(i-2:i+2,j) ! fine as long has nx_halo >= 2
4394
4395 flux_diff = 0
4396
4397 ! 1ST DO LEFT FACE
4398
4399 if (cs%u_face_mask(i-1,j) == 4.) then
4400
4401 flux_diff = flux_diff + g%dyCu(i-1,j) * time_step * cs%u_flux_bdry_val(i-1,j) * &
4402 cs%t_bdry_val(i-1,j) / g%areaT(i,j)
4403 else
4404
4405 ! get u-velocity at center of left face
4406 u_face = 0.5 * (cs%u_shelf(i-1,j-1) + cs%u_shelf(i-1,j))
4407
4408 if (u_face > 0) then !flux is into cell - we need info from h(i-2), h(i-1) if available
4409
4410 ! i may not cover all the cases.. but i cover the realistic ones
4411
4412 if (at_west_bdry .AND. (hmask(i-1,j) == 3)) then ! at western bdry but there is a
4413 ! thickness bdry condition, and the stencil contains it
4414 flux_diff = flux_diff + abs(u_face) * g%dyCu(i-1,j) * time_step * stencil(-1) / g%areaT(i,j)
4415
4416 elseif (hmask(i-1,j) * hmask(i-2,j) == 1) then ! h(i-2) and h(i-1) are valid
4417 phi = slope_limiter(stencil(-1)-stencil(-2), stencil(0)-stencil(-1))
4418 flux_diff = flux_diff + ((abs(u_face) * g%dyCu(i-1,j)* time_step / g%areaT(i,j)) * &
4419 (stencil(-1) - (phi * (stencil(-1)-stencil(0))/2)))
4420
4421 else ! h(i-1) is valid
4422 ! (o.w. flux would most likely be out of cell)
4423 ! but h(i-2) is not
4424
4425 flux_diff = flux_diff + abs(u_face) * g%dyCu(i-1,j) * time_step / g%areaT(i,j) * stencil(-1)
4426
4427 endif
4428
4429 elseif (u_face < 0) then !flux is out of cell - we need info from h(i-1), h(i+1) if available
4430 if (hmask(i-1,j) * hmask(i+1,j) == 1) then ! h(i-1) and h(i+1) are both valid
4431 phi = slope_limiter(stencil(0)-stencil(1), stencil(-1)-stencil(0))
4432 flux_diff = flux_diff - ((abs(u_face) * g%dyCu(i-1,j) * time_step / g%areaT(i,j)) * &
4433 (stencil(0) - (phi * (stencil(0)-stencil(-1))/2)))
4434
4435 else
4436 flux_diff = flux_diff - abs(u_face) * g%dyCu(i-1,j) * time_step / g%areaT(i,j) * stencil(0)
4437 endif
4438 endif
4439 endif
4440
4441 ! NEXT DO RIGHT FACE
4442
4443 ! get u-velocity at center of eastern face
4444
4445 if (cs%u_face_mask(i,j) == 4.) then
4446
4447 flux_diff = flux_diff + g%dyCu(i,j) * time_step * cs%u_flux_bdry_val(i,j) *&
4448 cs%t_bdry_val(i+1,j) / g%areaT(i,j)
4449 else
4450
4451 u_face = 0.5 * (cs%u_shelf(i,j-1) + cs%u_shelf(i,j))
4452
4453 if (u_face < 0) then !flux is into cell - we need info from h(i+2), h(i+1) if available
4454
4455 if (at_east_bdry .AND. (hmask(i+1,j) == 3)) then ! at eastern bdry but there is a
4456 ! thickness bdry condition, and the stencil contains it
4457
4458 flux_diff = flux_diff + abs(u_face) * g%dyCu(i,j) * time_step * stencil(1) / g%areaT(i,j)
4459
4460 elseif (hmask(i+1,j) * hmask(i+2,j) == 1) then ! h(i+2) and h(i+1) are valid
4461
4462 phi = slope_limiter(stencil(1)-stencil(2), stencil(0)-stencil(1))
4463 flux_diff = flux_diff + ((abs(u_face) * g%dyCu(i,j) * time_step / g%areaT(i,j)) * &
4464 (stencil(1) - (phi * (stencil(1)-stencil(0))/2)))
4465
4466 else ! h(i+1) is valid
4467 ! (o.w. flux would most likely be out of cell)
4468 ! but h(i+2) is not
4469
4470 flux_diff = flux_diff + abs(u_face) * g%dyCu(i,j) * time_step / g%areaT(i,j) * stencil(1)
4471
4472 endif
4473
4474 elseif (u_face > 0) then !flux is out of cell - we need info from h(i-1), h(i+1) if available
4475
4476 if (hmask(i-1,j) * hmask(i+1,j) == 1) then ! h(i-1) and h(i+1) are both valid
4477
4478 phi = slope_limiter(stencil(0)-stencil(-1), stencil(1)-stencil(0))
4479 flux_diff = flux_diff - ((abs(u_face) * g%dyCu(i,j) * time_step / g%areaT(i,j)) * &
4480 (stencil(0) - (phi * (stencil(0)-stencil(1))/2)))
4481
4482 else ! h(i+1) is valid (o.w. flux would most likely be out of cell) but h(i+2) is not
4483
4484 flux_diff = flux_diff - abs(u_face) * g%dyCu(i,j) * time_step / g%areaT(i,j) * stencil(0)
4485
4486 endif
4487
4488 endif
4489
4490 h_after_uflux(i,j) = h_after_uflux(i,j) + flux_diff
4491
4492 endif
4493
4494 endif
4495
4496 endif
4497
4498 enddo ! i loop
4499
4500 endif
4501
4502 enddo ! j loop
4503
4504end subroutine ice_shelf_advect_temp_x
4505
4506subroutine ice_shelf_advect_temp_y(CS, G, time_step, hmask, h_after_uflux, h_after_vflux)
4507 type(ice_shelf_dyn_cs), intent(in) :: CS !< A pointer to the ice shelf control structure
4508 type(ocean_grid_type), intent(in) :: G !< The grid structure used by the ice shelf.
4509 real, intent(in) :: time_step !< The time step for this update [T ~> s].
4510 real, dimension(SZDI_(G),SZDJ_(G)), &
4511 intent(in) :: hmask !< A mask indicating which tracer points are
4512 !! partly or fully covered by an ice-shelf
4513 real, dimension(SZDI_(G),SZDJ_(G)), &
4514 intent(in) :: h_after_uflux !< The ice shelf thicknesses times temperature after
4515 !! the zonal mass fluxes [C Z ~> degC m].
4516 real, dimension(SZDI_(G),SZDJ_(G)), &
4517 intent(inout) :: h_after_vflux !< The ice shelf thicknesses times temperature after
4518 !! the meridional mass fluxes [C Z ~> degC m]
4519
4520 ! use will be made of ISS%hmask here - its value at the boundary will be zero, just like uncovered cells
4521 ! if there is an input bdry condition, the thickness there will be set in initialization
4522
4523 integer :: i, j, is, ie, js, je, isd, ied, jsd, jed
4524 integer :: i_off, j_off
4525 logical :: at_north_bdry, at_south_bdry
4526 real, dimension(-2:2) :: stencil ! A copy of the neighboring thicknesses times temperatures [C Z ~> degC m]
4527 real :: v_face ! Pseudo-meridional velocity at a cell face, positive if out [L T-1 ~> m s-1]
4528 real :: flux_diff ! The difference in fluxes [C Z ~> degC m]
4529 real :: phi
4530
4531 is = g%isc ; ie = g%iec ; js = g%jsc-1 ; je = g%jec+1 ; isd = g%isd ; ied = g%ied ; jsd = g%jsd ; jed = g%jed
4532 i_off = g%idg_offset ; j_off = g%jdg_offset
4533
4534 do i=isd+2,ied-2
4535 if (((i+i_off) <= g%domain%niglobal+g%domain%nihalo) .AND. &
4536 ((i+i_off) >= g%domain%nihalo+1)) then ! based on mehmet's code - only if btw east & west boundaries
4537
4538 stencil(:) = 0.0 ! This is probably unnecessary, as the code is written
4539
4540 do j=js,je
4541
4542 if (((j+j_off) <= g%domain%njglobal+g%domain%njhalo) .AND. &
4543 ((j+j_off) >= g%domain%njhalo+1)) then
4544
4545 if (j+j_off == g%domain%njhalo+1) then
4546 at_south_bdry=.true.
4547 else
4548 at_south_bdry=.false.
4549 endif
4550 if (j+j_off == g%domain%njglobal+g%domain%njhalo) then
4551 at_north_bdry=.true.
4552 else
4553 at_north_bdry=.false.
4554 endif
4555
4556 if (hmask(i,j) == 1) then
4557 h_after_vflux(i,j) = h_after_uflux(i,j)
4558
4559 stencil(:) = h_after_uflux(i,j-2:j+2) ! fine as long has ny_halo >= 2
4560 flux_diff = 0
4561
4562 ! 1ST DO south FACE
4563
4564 if (cs%v_face_mask(i,j-1) == 4.) then
4565
4566 flux_diff = flux_diff + g%dxCv(i,j-1) * time_step * cs%v_flux_bdry_val(i,j-1) * &
4567 cs%t_bdry_val(i,j-1)/ g%areaT(i,j)
4568 else
4569
4570 ! get u-velocity at center of west face
4571 v_face = 0.5 * (cs%v_shelf(i-1,j-1) + cs%v_shelf(i,j-1))
4572
4573 if (v_face > 0) then !flux is into cell - we need info from h(j-2), h(j-1) if available
4574
4575 ! i may not cover all the cases.. but i cover the realistic ones
4576
4577 if (at_south_bdry .AND. (hmask(i,j-1) == 3)) then ! at western bdry but there is a
4578 ! thickness bdry condition, and the stencil contains it
4579 flux_diff = flux_diff + abs(v_face) * g%dxCv(i,j-1) * time_step * stencil(-1) / g%areaT(i,j)
4580
4581 elseif (hmask(i,j-1) * hmask(i,j-2) == 1) then ! h(j-2) and h(j-1) are valid
4582
4583 phi = slope_limiter(stencil(-1)-stencil(-2), stencil(0)-stencil(-1))
4584 flux_diff = flux_diff + ((abs(v_face) * g%dxCv(i,j-1) * time_step / g%areaT(i,j)) * &
4585 (stencil(-1) - (phi * (stencil(-1)-stencil(0))/2)))
4586
4587 else ! h(j-1) is valid
4588 ! (o.w. flux would most likely be out of cell)
4589 ! but h(j-2) is not
4590 flux_diff = flux_diff + abs(v_face) * g%dxCv(i,j-1) * time_step / g%areaT(i,j) * stencil(-1)
4591 endif
4592
4593 elseif (v_face < 0) then !flux is out of cell - we need info from h(j-1), h(j+1) if available
4594
4595 if (hmask(i,j-1) * hmask(i,j+1) == 1) then ! h(j-1) and h(j+1) are both valid
4596 phi = slope_limiter(stencil(0)-stencil(1), stencil(-1)-stencil(0))
4597 flux_diff = flux_diff - ((abs(v_face) * g%dxCv(i,j-1) * time_step / g%areaT(i,j)) * &
4598 (stencil(0) - (phi * (stencil(0)-stencil(-1))/2)))
4599 else
4600 flux_diff = flux_diff - abs(v_face) * g%dxCv(i,j-1) * time_step / g%areaT(i,j) * stencil(0)
4601 endif
4602
4603 endif
4604
4605 endif
4606
4607 ! NEXT DO north FACE
4608
4609 if (cs%v_face_mask(i,j) == 4.) then
4610 flux_diff = flux_diff + g%dxCv(i,j) * time_step * cs%v_flux_bdry_val(i,j) *&
4611 cs%t_bdry_val(i,j+1)/ g%areaT(i,j)
4612 else
4613
4614 ! get u-velocity at center of east face
4615 v_face = 0.5 * (cs%v_shelf(i-1,j) + cs%v_shelf(i,j))
4616
4617 if (v_face < 0) then !flux is into cell - we need info from h(j+2), h(j+1) if available
4618
4619 if (at_north_bdry .AND. (hmask(i,j+1) == 3)) then ! at eastern bdry but there is a
4620 ! thickness bdry condition, and the stencil contains it
4621 flux_diff = flux_diff + abs(v_face) * g%dxCv(i,j) * time_step * stencil(1) / g%areaT(i,j)
4622 elseif (hmask(i,j+1) * hmask(i,j+2) == 1) then ! h(j+2) and h(j+1) are valid
4623 phi = slope_limiter(stencil(1)-stencil(2), stencil(0)-stencil(1))
4624 flux_diff = flux_diff + ((abs(v_face) * g%dxCv(i,j) * time_step / g%areaT(i,j)) * &
4625 (stencil(1) - (phi * (stencil(1)-stencil(0))/2)))
4626 else ! h(j+1) is valid
4627 ! (o.w. flux would most likely be out of cell)
4628 ! but h(j+2) is not
4629 flux_diff = flux_diff + abs(v_face) * g%dxCv(i,j) * time_step / g%areaT(i,j) * stencil(1)
4630 endif
4631
4632 elseif (v_face > 0) then !flux is out of cell - we need info from h(j-1), h(j+1) if available
4633
4634 if (hmask(i,j-1) * hmask(i,j+1) == 1) then ! h(j-1) and h(j+1) are both valid
4635 phi = slope_limiter(stencil(0)-stencil(-1), stencil(1)-stencil(0))
4636 flux_diff = flux_diff - ((abs(v_face) * g%dxCv(i,j) * time_step / g%areaT(i,j)) * &
4637 (stencil(0) - (phi * (stencil(0)-stencil(1))/2)))
4638 else ! h(j+1) is valid
4639 ! (o.w. flux would most likely be out of cell)
4640 ! but h(j+2) is not
4641 flux_diff = flux_diff - abs(v_face) * g%dxCv(i,j) * time_step / g%areaT(i,j) * stencil(0)
4642 endif
4643
4644 endif
4645
4646 endif
4647
4648 h_after_vflux(i,j) = h_after_vflux(i,j) + flux_diff
4649 endif
4650 endif
4651 enddo ! j loop
4652 endif
4653 enddo ! i loop
4654
4655end subroutine ice_shelf_advect_temp_y
4656
4657end module mom_ice_shelf_dynamics