Skip to content

Commit defa661

Browse files
committed
Add test for composite rheology
1 parent 9047aa2 commit defa661

File tree

1 file changed

+377
-0
lines changed

1 file changed

+377
-0
lines changed

test/composite_rheology_test.py

Lines changed: 377 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,377 @@
1+
import pytest
2+
import numpy as np
3+
import sympy
4+
from sympy import lambdify, simplify, diff
5+
import firedrake
6+
from firedrake import (
7+
as_vector,
8+
max_value,
9+
Constant,
10+
Function,
11+
derivative,
12+
NonlinearVariationalProblem,
13+
NonlinearVariationalSolver,
14+
DirichletBC,
15+
)
16+
from icepack2 import constants
17+
from icepack2.model.minimization import (
18+
viscous_power,
19+
friction_power,
20+
calving_terminus,
21+
momentum_balance,
22+
ice_shelf_momentum_balance,
23+
)
24+
25+
26+
sparams = {
27+
"solver_parameters": {
28+
"snes_type": "newtonls",
29+
"snes_max_it": 200,
30+
"snes_linesearch_type": "nleqerr",
31+
"ksp_type": "gmres",
32+
"pc_type": "lu",
33+
"pc_factor_mat_solver_type": "mumps",
34+
},
35+
}
36+
37+
38+
@pytest.mark.parametrize("degree", (1, 2))
39+
def test_composite_rheology_floating(degree):
40+
lx, ly = 20e3, 20e3
41+
Lx, Ly = Constant(lx), Constant(ly)
42+
h0, dh = Constant(500.0), Constant(100.0)
43+
u_inflow = Constant(100.0)
44+
45+
ρ_I = Constant(constants.ice_density)
46+
ρ_W = Constant(constants.water_density)
47+
g = Constant(constants.gravity)
48+
49+
# Use a combination of two different rheology exponents.
50+
n_1 = firedrake.Constant(1.0)
51+
n_3 = firedrake.Constant(constants.glen_flow_law)
52+
53+
τ_c = Constant(0.1)
54+
ε_c = Constant(0.01)
55+
A_1 = ε_c / τ_c ** n_1
56+
A_3 = ε_c / τ_c ** n_3
57+
58+
# This is an exact solution for the velocity of a floating ice shelf with
59+
# constant temperature and linearly decreasing thickness. See Greve and
60+
# Blatter for the derivation.
61+
def exact_δu(x, n, A):
62+
ρ = ρ_I * (1 - ρ_I / ρ_W)
63+
h = h0 - dh * x / Lx
64+
P = ρ * g * h / 4
65+
P_0 = ρ * g * h0 / 4
66+
δP = ρ * g * dh / 4
67+
return Lx * A * (P_0 ** (n + 1) - P ** (n + 1)) / ((n + 1) * δP)
68+
69+
errors, mesh_sizes = [], []
70+
k_min, k_max, num_steps = 5 - degree, 8 - degree, 9
71+
for nx in np.logspace(k_min, k_max, num_steps, base=2, dtype=int):
72+
mesh = firedrake.RectangleMesh(nx, nx, lx, ly, diagonal="crossed")
73+
x, y = firedrake.SpatialCoordinate(mesh)
74+
75+
# Make some function spaces.
76+
cg = firedrake.FiniteElement("CG", "triangle", degree)
77+
dg = firedrake.FiniteElement("DG", "triangle", degree - 1)
78+
Q = firedrake.FunctionSpace(mesh, cg)
79+
V = firedrake.VectorFunctionSpace(mesh, cg)
80+
Σ = firedrake.TensorFunctionSpace(mesh, dg, symmetry=True)
81+
Z = V * Σ
82+
z = Function(Z)
83+
z.sub(0).assign(Constant(u_inflow, 0))
84+
85+
expr = u_inflow + exact_δu(x, n_1, A_1) + exact_δu(x, n_3, A_3)
86+
u_exact = Function(V).interpolate(as_vector((expr, 0)))
87+
88+
h = Function(Q).interpolate(h0 - dh * x / Lx)
89+
90+
# Make the boundary conditions.
91+
inflow_ids = (1,)
92+
outflow_ids = (2,)
93+
side_wall_ids = (3, 4)
94+
95+
inflow_bc = DirichletBC(Z.sub(0), Constant((u_inflow, 0)), inflow_ids)
96+
side_wall_bc = DirichletBC(Z.sub(0).sub(1), 0, side_wall_ids)
97+
bcs = [inflow_bc, side_wall_bc]
98+
99+
# Make the model specification and solver.
100+
u, M = firedrake.split(z)
101+
fields = {"velocity": u, "membrane_stress": M, "thickness": h}
102+
boundary_ids = {"outflow_ids": outflow_ids}
103+
104+
L = (
105+
viscous_power(**fields, flow_law_exponent=n_1, flow_law_coefficient=A_1) +
106+
viscous_power(**fields, flow_law_exponent=n_3, flow_law_coefficient=A_3) +
107+
ice_shelf_momentum_balance(**fields, **boundary_ids)
108+
)
109+
110+
F = derivative(L, z)
111+
qdegree = max(8, degree ** constants.glen_flow_law)
112+
pparams = {"form_compiler_parameters": {"quadrature_degree": qdegree}}
113+
problem = NonlinearVariationalProblem(F, z, bcs, **pparams)
114+
solver = NonlinearVariationalSolver(problem, **sparams)
115+
116+
# Solve the problem using a continuation method -- step the flow law
117+
# exponent from 1 to 3 in small intervals.
118+
num_steps = 5
119+
for exponent in np.linspace(1.0, constants.glen_flow_law, num_steps):
120+
n_3.assign(exponent)
121+
solver.solve()
122+
123+
u, M = z.subfunctions
124+
error = firedrake.norm(u - u_exact) / firedrake.norm(u_exact)
125+
δx = mesh.cell_sizes.dat.data_ro.min()
126+
mesh_sizes.append(δx)
127+
errors.append(error)
128+
print(".", end="", flush=True)
129+
130+
log_mesh_sizes = np.log2(np.array(mesh_sizes))
131+
log_errors = np.log2(np.array(errors))
132+
slope, intercept = np.polyfit(log_mesh_sizes, log_errors, 1)
133+
print(f"degree {degree}: log(error) ~= {slope:g} * log(dx) {intercept:+g}")
134+
assert slope > degree + 0.9
135+
136+
137+
def test_manufactured_solution():
138+
L, ρ_I, ρ_W, g = sympy.symbols("L ρ_I ρ_W g", real=True, positive=True)
139+
A_1, A_3 = sympy.symbols("A_1 A_3", real=True, positive=True)
140+
n_1, n_3 = sympy.symbols("n_1 n_3", integer=True, positive=True)
141+
142+
def strain_rate(M, A_1, A_3, n_1, n_3):
143+
return A_1 * (M / 2)**n_1 + A_3 * (M / 2)**n_3
144+
145+
def sliding_velocity(τ, K_1, K_3, m_1, m_3):
146+
return -K_1 * τ**m_1 - K_3 * τ**m_3
147+
148+
def stress_balance(x, h, s, M, τ):
149+
return diff(h * M, x) + τ - ρ_I * g * h * diff(s, x)
150+
151+
def boundary_condition(x, u, h, s, M):
152+
d = (s - h).subs(x, L)
153+
τ_c = (ρ_I * g * h**2 - ρ_W * g * d**2) / 2
154+
return (h * M - τ_c).subs(x, L)
155+
156+
# Create the thickness and surface elevation
157+
x = sympy.symbols("x", real=True)
158+
h0, δh = sympy.symbols("h0 δh", real=True, positive=True)
159+
h = h0 - δh * x / L
160+
161+
s0, δs = sympy.symbols("s0 δs", real=True, positive=True)
162+
s = s0 - δs * x / L
163+
164+
# The variable `β` is the fraction of the driving stress that membrane
165+
# stress divergence will take up; this is chosen to eliminate some extra
166+
# terms and make the rest of the algebra simpler. See the website also.
167+
h_L = h0 - δh
168+
s_L = s0 - δs
169+
β = δh / δs * (ρ_I * h_L ** 2 - ρ_W * (s_L - h_L) ** 2) / (ρ_I * h_L**2)
170+
171+
ρ = β * ρ_I * δs / δh
172+
P = ρ * g * h / 4
173+
δP = ρ * g * δh / 4
174+
P0 = ρ * g * h0 / 4
175+
176+
M = ρ * g * h / 2
177+
τ = -(1 - β) * ρ_I * g * h * δs / L
178+
179+
# The exact velocity field looks similar to the floating case by using a
180+
# careful choice of `β`.
181+
u_0 = sympy.symbols("u_0", real=True, positive=True)
182+
δu_1 = L * A_1 * (P0 ** (n_1 + 1) - P ** (n_1 + 1)) / ((n_1 + 1) * δP)
183+
δu_3 = L * A_3 * (P0 ** (n_3 + 1) - P ** (n_3 + 1)) / ((n_3 + 1) * δP)
184+
u = u_0 + δu_1 + δu_3
185+
186+
# Make the linear part of the rheology a much smaller effect than the
187+
# nonlinear part -- at a stress of 100kPa, the nonlinear rheology gives
188+
# a strain rate of 10 m/yr/km, while the linear part only 1 m/yr/km.
189+
ε_1 = 0.001
190+
ε_3 = 0.01
191+
τ_c = 0.1
192+
193+
values = {
194+
u_0: 100,
195+
h0: 500,
196+
δh: 100,
197+
s0: 150,
198+
δs: 90,
199+
L: 20e3,
200+
A_1: ε_1 / τ_c,
201+
A_3: ε_3 / τ_c ** constants.glen_flow_law,
202+
ρ_I: constants.ice_density,
203+
ρ_W: constants.water_density,
204+
n_1: 1,
205+
n_3: constants.glen_flow_law,
206+
g: constants.gravity,
207+
}
208+
209+
# Check the momentum conservation law and boundary condition.
210+
τ_b = lambdify(x, τ.subs(values), "numpy")
211+
τ_d = lambdify(x, (-ρ_I * g * h * diff(s, x)).subs(values), "numpy")
212+
τ_m = lambdify(x, simplify(diff(h * M, x)).subs(values), "numpy")
213+
xs = np.linspace(0, values[L], 21)
214+
215+
tolerance = 1e-8
216+
assert abs(boundary_condition(x, u, h, s, M).subs(values)) < tolerance
217+
218+
stress_norm = np.max(np.abs(τ_d(xs)))
219+
assert np.max(np.abs(τ_m(xs) + τ_b(xs) + τ_d(xs))) < tolerance * stress_norm
220+
221+
# Check the constitutive relation.
222+
ε_u = lambdify(x, simplify(diff(u, x).subs(values)), "numpy")
223+
ε_M = lambdify(
224+
x, simplify(strain_rate(M, A_1, A_3, n_1, n_3).subs(values)), "numpy"
225+
)
226+
assert np.max(np.abs(ε_u(xs) - ε_M(xs))) < tolerance * np.max(np.abs(ε_u(xs)))
227+
228+
# Check the sliding law.
229+
α = 0.1
230+
m_1 = 1
231+
m_3 = constants.glen_flow_law
232+
K_1 = α * u / abs(τ) ** m_1
233+
K_3 = (u - K_1 * abs(τ) ** m_1) / abs(τ)**m_3
234+
235+
U = lambdify(x, simplify(u).subs(values), "numpy")
236+
U_b = lambdify(
237+
x, simplify(sliding_velocity(τ, K_1, K_3, m_1, m_3).subs(values)), "numpy"
238+
)
239+
assert np.max(np.abs(U(xs) - U_b(xs))) < tolerance * np.max(np.abs(U(xs)))
240+
241+
242+
@pytest.mark.parametrize("degree", (1, 2))
243+
def test_composite_rheology_grounded(degree):
244+
lx, ly = 20e3, 20e3
245+
Lx, Ly = Constant(lx), Constant(ly)
246+
h0, dh = Constant(500.0), Constant(100.0)
247+
s0, ds = Constant(150.0), Constant(90.0)
248+
u_inflow = Constant(100.0)
249+
250+
# See the previous test for explanations of all the setup.
251+
ρ_I = Constant(constants.ice_density)
252+
ρ_W = Constant(constants.water_density)
253+
g = Constant(constants.gravity)
254+
255+
n_1 = firedrake.Constant(1.0)
256+
n_3 = firedrake.Constant(constants.glen_flow_law)
257+
258+
m_1 = firedrake.Constant(1.0)
259+
m_3 = firedrake.Constant(constants.weertman_sliding_law)
260+
261+
τ_c = Constant(0.1)
262+
ε_1 = Constant(0.001)
263+
ε_3 = Constant(0.01)
264+
A_1 = ε_1 / τ_c ** n_1
265+
A_3 = ε_3 / τ_c ** n_3
266+
267+
h_L = h0 - dh
268+
s_L = s0 - ds
269+
β = dh / ds * (ρ_I * h_L**2 - ρ_W * (s_L - h_L)**2) / (ρ_I * h_L**2)
270+
271+
def exact_δu(x, n, A):
272+
ρ = β * ρ_I * ds / dh
273+
h = h0 - dh * x / Lx
274+
P = ρ * g * h / 4
275+
dP = ρ * g * dh / 4
276+
P0 = ρ * g * h0 / 4
277+
du = Lx * A * (P0 ** (n + 1) - P ** (n + 1)) / ((n + 1) * dP)
278+
return du
279+
280+
errors, mesh_sizes = [], []
281+
k_min, k_max, num_steps = 5 - degree, 8 - degree, 9
282+
for nx in np.logspace(k_min, k_max, num_steps, base=2, dtype=int):
283+
mesh = firedrake.RectangleMesh(nx, nx, lx, ly, diagonal="crossed")
284+
x, y = firedrake.SpatialCoordinate(mesh)
285+
286+
# Make some function spaces.
287+
cg = firedrake.FiniteElement("CG", "triangle", degree)
288+
dg = firedrake.FiniteElement("DG", "triangle", degree - 1)
289+
Q = firedrake.FunctionSpace(mesh, cg)
290+
V = firedrake.VectorFunctionSpace(mesh, cg)
291+
Σ = firedrake.TensorFunctionSpace(mesh, dg, symmetry=True)
292+
T = firedrake.VectorFunctionSpace(mesh, dg)
293+
Z = V * Σ * T
294+
z = Function(Z)
295+
z.sub(0).assign(Constant((u_inflow, 0)))
296+
297+
# Make the exact velocity, thickness, surface elevation, and sliding
298+
# coefficients.
299+
u_expr = u_inflow + exact_δu(x, n_1, A_1) + exact_δu(x, n_3, A_3)
300+
u_exact = Function(V).interpolate(as_vector((u_expr, 0)))
301+
302+
h = Function(Q).interpolate(h0 - dh * x / Lx)
303+
s = Function(Q).interpolate(s0 - ds * x / Lx)
304+
305+
α = Constant(0.1)
306+
τ_expr = (1 - β) * ρ_I * g * h * ds / Lx
307+
K_1 = α * u_expr / τ_expr ** m_1
308+
K_3 = (u_expr - K_1 * τ_expr ** m_1) / τ_expr ** m_3
309+
310+
# Create the boundary conditions.
311+
inflow_ids = (1,)
312+
outflow_ids = (2,)
313+
side_wall_ids = (3, 4)
314+
315+
inflow_bc = DirichletBC(Z.sub(0), Constant((u_inflow, 0)), inflow_ids)
316+
side_wall_bc = DirichletBC(Z.sub(0).sub(1), 0, side_wall_ids)
317+
bcs = [inflow_bc, side_wall_bc]
318+
319+
# Create the model specification and solvers.
320+
u, M, τ = firedrake.split(z)
321+
fields = {
322+
"velocity": u,
323+
"membrane_stress": M,
324+
"basal_stress": τ,
325+
"thickness": h,
326+
"surface": s,
327+
}
328+
329+
rheology_1 = {
330+
"flow_law_exponent": n_1,
331+
"flow_law_coefficient": A_1,
332+
"sliding_exponent": m_1,
333+
"sliding_coefficient": K_1,
334+
}
335+
336+
rheology_3 = {
337+
"flow_law_exponent": n_3,
338+
"flow_law_coefficient": A_3,
339+
"sliding_exponent": m_3,
340+
"sliding_coefficient": K_3,
341+
}
342+
343+
boundary_ids = {"outflow_ids": outflow_ids}
344+
L = (
345+
viscous_power(**fields, **rheology_1) +
346+
viscous_power(**fields, **rheology_3) +
347+
friction_power(**fields, **rheology_1) +
348+
friction_power(**fields, **rheology_3) +
349+
calving_terminus(**fields, **boundary_ids) +
350+
momentum_balance(**fields)
351+
)
352+
353+
F = derivative(L, z)
354+
qdegree = max(8, degree ** constants.glen_flow_law)
355+
pparams = {"form_compiler_parameters": {"quadrature_degree": qdegree}}
356+
problem = NonlinearVariationalProblem(F, z, bcs, **pparams)
357+
solver = NonlinearVariationalSolver(problem, **sparams)
358+
359+
# Solve the problem using a continuation method.
360+
num_steps = 5
361+
for exponent in np.linspace(1.0, constants.glen_flow_law, num_steps):
362+
n_3.assign(exponent)
363+
m_3.assign(exponent)
364+
solver.solve()
365+
366+
u, M, τ = z.subfunctions
367+
error = firedrake.norm(u - u_exact) / firedrake.norm(u_exact)
368+
δx = mesh.cell_sizes.dat.data_ro.min()
369+
mesh_sizes.append(δx)
370+
errors.append(error)
371+
print(".", end="", flush=True)
372+
373+
log_mesh_sizes = np.log2(np.array(mesh_sizes))
374+
log_errors = np.log2(np.array(errors))
375+
slope, intercept = np.polyfit(log_mesh_sizes, log_errors, 1)
376+
print(f"degree {degree}: log(error) ~= {slope:g} * log(dx) {intercept:+g}")
377+
assert slope > degree + 0.9

0 commit comments

Comments
 (0)