def side_friction(**kwargs): r"""Return the side wall friction part of the action functional The component of the action functional due to friction along the side walls of the domain is .. math:: E(u) = -\frac{m}{m + 1}\int_\Gamma h\tau(u, C_s)\cdot u\; ds where :math:`\tau(u, C_s)` is the side wall shear stress, :math:`ds` is the element of surface area and :math:`\Gamma` are the side walls. Side wall friction is relevant for glaciers that flow through a fjord with rock walls on either side. """ u, h = get_kwargs_alt(kwargs, ('velocity', 'thickness'), ('u', 'h')) Cs = kwargs.get('side_friction', kwargs.get('Cs', firedrake.Constant(0.))) mesh = u.ufl_domain() if mesh.geometric_dimension() == 2: ν = firedrake.FacetNormal(mesh) else: ν = facet_normal_2(mesh) u_t = u - inner(u, ν) * ν τ = friction_stress(u_t, Cs) return -m / (m + 1) * h * inner(τ, u_t)
def normal_flow_penalty(**kwargs): r"""Return the penalty for flow normal to the domain boundary For problems where a glacier flows along some boundary, e.g. a fjord wall, the velocity has to be parallel to this boundary. Rather than enforce this boundary condition directly, we add a penalty for normal flow to the action functional. """ u, = get_kwargs_alt(kwargs, ('velocity', ), ('u', )) scale = kwargs.get('scale', firedrake.Constant(1.)) mesh = u.ufl_domain() if mesh.geometric_dimension() == 2: ν = firedrake.FacetNormal(mesh) elif mesh.geometric_dimension() == 3: ν = facet_normal_2(mesh) L = diameter(mesh) δx = firedrake.FacetArea(mesh) # Get the polynomial degree in the horizontal direction of the velocity # field -- if extruded, the element degree is a tuple of the horizontal # and vertical degrees. degree = u.ufl_function_space().ufl_element().degree() if isinstance(degree, tuple): d = degree[0] else: d = degree exponent = kwargs.get('exponent', d + 1) penalty = scale * (L / δx)**exponent return 0.5 * penalty * inner(u, ν)**2
def advective_flux(self, **kwargs): keys = ('energy', 'velocity', 'vertical_velocity', 'thickness', 'energy_inflow', 'energy_surface') keys_alt = ('E', 'u', 'w', 'h', 'E_inflow', 'E_surface') E, u, w, h, E_inflow, E_surface = get_kwargs_alt( kwargs, keys, keys_alt) Q = E.function_space() ψ = firedrake.TestFunction(Q) U = firedrake.as_vector((u[0], u[1], w)) flux_cells = -E * inner(U, grad(ψ)) * h * dx mesh = Q.mesh() ν = facet_normal_2(mesh) outflow = firedrake.max_value(inner(u, ν), 0) inflow = firedrake.min_value(inner(u, ν), 0) flux_outflow = (E * outflow * ψ * h * ds_v + E * firedrake.max_value(-w, 0) * ψ * h * ds_b + E * firedrake.max_value(+w, 0) * ψ * h * ds_t) flux_inflow = (E_inflow * inflow * ψ * h * ds_v + E_surface * firedrake.min_value(-w, 0) * ψ * h * ds_b + E_surface * firedrake.min_value(+w, 0) * ψ * h * ds_t) return flux_cells + flux_outflow + flux_inflow
def terminus(u, h, s): r"""Return the terminal stress part of the hybrid model action functional The power exerted due to stress at the calving terminus :math:`\Gamma` is .. math:: E(u) = \int_\Gamma\int_0^1\left(\rho_Ig(1 - \zeta) - \rho_Wg(\zeta_{\text{sl}} - \zeta)_+\right)u\cdot\nu\; h\, d\zeta\; ds where :math:`\zeta_\text{sl}` is the relative depth to sea level and the :math:`(\zeta_\text{sl} - \zeta)_+` denotes only the positive part. Parameters ---------- u : firedrake.Function ice velocity h : firedrake.Function ice thickness s : firedrake.Function ice surface elevation ice_front_ids : list of int numeric IDs of the parts of the boundary corresponding to the calving front """ mesh = u.ufl_domain() zdegree = u.ufl_element().degree()[1] x, y, ζ = firedrake.SpatialCoordinate(mesh) b = s - h ζ_sl = firedrake.max_value(-b, 0) / h p_W = ρ_W * g * h * _pressure_approx(zdegree + 1)(ζ, ζ_sl) p_I = ρ_I * g * h * (1 - ζ) ν = facet_normal_2(mesh) return (p_I - p_W) * inner(u, ν) * h
def _advect(self, dt, E, u, w, h, s, E_inflow, E_surface): Q = E.function_space() φ, ψ = firedrake.TrialFunction(Q), firedrake.TestFunction(Q) U = firedrake.as_vector((u[0], u[1], w)) flux_cells = -φ * inner(U, grad(ψ)) * h * dx mesh = Q.mesh() ν = facet_normal_2(mesh) outflow = firedrake.max_value(inner(u, ν), 0) inflow = firedrake.min_value(inner(u, ν), 0) flux_outflow = φ * ψ * outflow * h * ds_v + \ φ * ψ * firedrake.max_value(-w, 0) * h * ds_b + \ φ * ψ * firedrake.max_value(+w, 0) * h * ds_t F = φ * ψ * h * dx + dt * (flux_cells + flux_outflow) flux_inflow = -E_inflow * ψ * inflow * h * ds_v \ -E_surface * ψ * firedrake.min_value(-w, 0) * h * ds_b \ -E_surface * ψ * firedrake.min_value(+w, 0) * h * ds_t A = E * ψ * h * dx + dt * flux_inflow solver_parameters = {'ksp_type': 'preonly', 'pc_type': 'lu'} degree_E = E.ufl_element().degree() degree_u = u.ufl_element().degree() degree = (3 * degree_E[0] + degree_u[0], 2 * degree_E[1] + degree_u[1]) form_compiler_parameters = {'quadrature_degree': degree} firedrake.solve(F == A, E, solver_parameters=solver_parameters, form_compiler_parameters=form_compiler_parameters)
def normal_flow_penalty(u, h, scale=1.0, exponent=None, side_wall_ids=()): r"""Return the penalty for flow normal to the domain boundary For problems where a glacier flows along some boundary, e.g. a fjord wall, the velocity has to be parallel to this boundary. Rather than enforce this boundary condition directly, we add a penalty for normal flow to the action functional. """ mesh = u.ufl_domain() ν = facet_normal_2(mesh) # Note that this quantity has units of [length] x [dimensionless] because # the mesh has a "thickness" of 1! If it had dimensions of physical # thickness, we would instead use the square root of the facet area. δx = firedrake.FacetArea(mesh) L = diameter(mesh) d = u.ufl_function_space().ufl_element().degree()[0] exponent = d + 1 if (exponent is None) else exponent penalty = scale * (L / δx)**exponent return 0.5 * penalty * inner(u, ν)**2 * h * ds_v(tuple(side_wall_ids))
def side_friction(u, h, Cs=firedrake.Constant(0), side_wall_ids=()): r"""Return the side wall friction part of the action functional The component of the action functional due to friction along the side walls of the domain is .. math:: E(u) = -\frac{m}{m + 1}\int_\Gamma\int_0^1 \tau(u, C_s)\cdot u \hspace{2pt}hd\zeta\hspace{2pt}ds where :math:`\tau(u, C_s)` is the side wall shear stress, :math:`ds` is the element of surface area and :math:`\Gamma` are the side walls. Side wall friction is relevant for glaciers that flow through a fjord with rock walls on either side. """ mesh = u.ufl_domain() ν = facet_normal_2(mesh) u_t = u - inner(u, ν) * ν τ = friction_stress(u_t, Cs) ids = tuple(side_wall_ids) return -m/(m + 1) * inner(τ, u_t) * h * ds_v(domain=mesh, subdomain_id=ids)
def terminus(u, h, s, ice_front_ids=()): r"""Return the terminal stress part of the hybrid model action functional The power exerted due to stress at the calving terminus :math:`\Gamma` is .. math:: E(u) = \int_\Gamma\int_0^1\left(\rho_Ig(1 - \zeta) - \rho_Wg(\zeta_{\text{sl}} - \zeta)_+\right)hd\zeta ds where :math:`\zeta_\text{sl}` is the relative depth to sea level and the :math:`(\zeta - \zeta_\text{sl})` denotes only the positive part. Parameters ---------- u : firedrake.Function ice velocity h : firedrake.Function ice thickness s : firedrake.Function ice surface elevation ice_front_ids : list of int numeric IDs of the parts of the boundary corresponding to the calving front """ xdegree_u, zdegree_u = u.ufl_element().degree() degree_h = h.ufl_element().degree()[0] degree = (xdegree_u + degree_h, 2 * zdegree_u + 1) metadata = {'quadrature_degree': degree} x, y, ζ = firedrake.SpatialCoordinate(u.ufl_domain()) b = s - h ζ_sl = firedrake.max_value(-b, 0) / h p_W = ρ_W * g * h * _pressure_approx(zdegree_u + 1)(ζ, ζ_sl) p_I = ρ_I * g * h * (1 - ζ) ν = facet_normal_2(u.ufl_domain()) dγ = ds_v(tuple(ice_front_ids), metadata=metadata) return (p_I - p_W) * inner(u, ν) * h * dγ
def solve(self, dt, h0, a, u, h_inflow=None): h_inflow = h_inflow if h_inflow is not None else h0 Q = h0.function_space() h, φ = firedrake.TrialFunction(Q), firedrake.TestFunction(Q) ν = facet_normal_2(Q.mesh()) outflow = firedrake.max_value(inner(u, ν), 0) inflow = firedrake.min_value(inner(u, ν), 0) flux_cells = -h * inner(u, grad_2(φ)) * dx flux_out = h * φ * outflow * ds_v F = h * φ * dx + dt * (flux_cells + flux_out) accumulation = a * φ * dx flux_in = -h_inflow * φ * inflow * ds_v A = h0 * φ * dx + dt * (accumulation + flux_in) h = h0.copy(deepcopy=True) solver_parameters = {'ksp_type': 'preonly', 'pc_type': 'lu'} firedrake.solve(F == A, h, solver_parameters=solver_parameters) return h