Skip to content

Commit 3a3a006

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

File tree

1 file changed

+376
-0
lines changed

1 file changed

+376
-0
lines changed

test/composite_rheology_test.py

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

0 commit comments

Comments
 (0)