-
Notifications
You must be signed in to change notification settings - Fork 13
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Tests for Terminator Toy and SplitPrescribedTransport #465
Merged
Merged
Changes from all commits
Commits
Show all changes
11 commits
Select commit
Hold shift + click to select a range
e433e77
SplitPrescribedTransport test of transporting velocity
ta440 730b107
preciptation integration test now uses SplitPrescribedTransport
ta440 4a2df00
added a test for the terminator physics
ta440 3076528
lint terminator physics test
ta440 2918b52
lint
ta440 b0b7750
Merge remote-tracking branch 'origin/main' into terminator_tests
ta440 78a0c0f
ref_level=2 for terminator test
ta440 216c7d9
lint
ta440 7275a70
lint
ta440 6200cc0
undo file deletion
ta440 366c7a5
Merge branch 'main' into terminator_tests
tommbendall File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
@@ -0,0 +1,110 @@ | ||
""" | ||
This tests the terminator toy physics scheme that models the interaction | ||
of two chemical species through coupled ODEs. | ||
""" | ||
|
||
from gusto import * | ||
from firedrake import IcosahedralSphereMesh, Constant, cos, \ | ||
sin, SpatialCoordinate, Function, max_value, as_vector, \ | ||
errornorm, norm | ||
import numpy as np | ||
|
||
|
||
def run_terminator_toy(dirname): | ||
|
||
# ------------------------------------------------------------------------ # | ||
# Set up model objects | ||
# ------------------------------------------------------------------------ # | ||
|
||
# A much larger timestep than in proper simulations, as this | ||
# tests moving towards a steady state with no flow. | ||
dt = 50000. | ||
|
||
# Make the mesh and domain | ||
R = 6371220. | ||
mesh = IcosahedralSphereMesh(radius=R, | ||
refinement_level=2, degree=2) | ||
|
||
# get lat lon coordinates | ||
x = SpatialCoordinate(mesh) | ||
lamda, theta, _ = lonlatr_from_xyz(x[0], x[1], x[2]) | ||
|
||
domain = Domain(mesh, dt, 'BDM', 1) | ||
|
||
# Define the interacting species | ||
X = ActiveTracer(name='X', space='DG', | ||
variable_type=TracerVariableType.mixing_ratio, | ||
transport_eqn=TransportEquationType.advective) | ||
|
||
X2 = ActiveTracer(name='X2', space='DG', | ||
variable_type=TracerVariableType.mixing_ratio, | ||
transport_eqn=TransportEquationType.advective) | ||
|
||
tracers = [X, X2] | ||
|
||
# Equation | ||
V = domain.spaces("HDiv") | ||
eqn = CoupledTransportEquation(domain, active_tracers=tracers, Vu=V) | ||
|
||
output = OutputParameters(dirname=dirname+"/terminator_toy", | ||
dumpfreq=10) | ||
io = IO(domain, output) | ||
|
||
# Define the reaction rates: | ||
theta_c = np.pi/9. | ||
lamda_c = -np.pi/3. | ||
|
||
k1 = max_value(0, sin(theta)*sin(theta_c) + cos(theta)*cos(theta_c)*cos(lamda-lamda_c)) | ||
k2 = 1 | ||
|
||
physics_schemes = [(TerminatorToy(eqn, k1=k1, k2=k2, species1_name='X', | ||
species2_name='X2'), BackwardEuler(domain))] | ||
|
||
# Set up a non-divergent, time-varying, velocity field | ||
def u_t(t): | ||
return as_vector([Constant(0)*lamda, Constant(0)*lamda, Constant(0)*lamda]) | ||
|
||
X_T_0 = 4e-6 | ||
X_0 = X_T_0 + 0*lamda | ||
X2_0 = 0*lamda | ||
|
||
transport_scheme = SSPRK3(domain) | ||
transport_method = [DGUpwind(eqn, 'X'), DGUpwind(eqn, 'X2')] | ||
|
||
stepper = SplitPrescribedTransport(eqn, transport_scheme, io, | ||
spatial_methods=transport_method, | ||
physics_schemes=physics_schemes, | ||
prescribed_transporting_velocity=u_t) | ||
|
||
stepper.fields("X").interpolate(X_0) | ||
stepper.fields("X2").interpolate(X2_0) | ||
|
||
stepper.run(t=0, tmax=10*dt) | ||
|
||
# Compute the steady state solution to compare to | ||
steady_space = domain.spaces('DG') | ||
X_steady = Function(steady_space) | ||
X2_steady = Function(steady_space) | ||
|
||
r = k1/(4*k2) | ||
D_val = sqrt(r**2 + 2*X_T_0*r) | ||
|
||
X_steady.interpolate(D_val - r) | ||
X2_steady.interpolate(0.5*(X_T_0 - D_val + r)) | ||
|
||
return stepper, X_steady, X2_steady | ||
|
||
|
||
def test_terminator_toy_setup(tmpdir): | ||
dirname = str(tmpdir) | ||
stepper, X_steady, X2_steady = run_terminator_toy(dirname) | ||
X_field = stepper.fields("X") | ||
X2_field = stepper.fields("X2") | ||
|
||
print(errornorm(X_field, X_steady)/norm(X_steady)) | ||
print(errornorm(X2_field, X2_steady)/norm(X2_steady)) | ||
|
||
# Assert that the physics scheme has sufficiently moved | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. This seems like a sensible test to me! |
||
# the species fields near their steady state solutions | ||
assert errornorm(X_field, X_steady)/norm(X_steady) < 0.4, "The X field is not sufficiently close to the steady state profile" | ||
assert errornorm(X2_field, X2_steady)/norm(X2_steady) < 0.4, "The X2 field is not sufficiently close to the steady state profile" |
This file contains bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
Original file line number | Diff line number | Diff line change |
---|---|---|
|
@@ -37,7 +37,7 @@ def run_wind_drag(dirname, implicit_formulation): | |
eqn = CompressibleEulerEquations(domain, parameters) | ||
|
||
# I/O | ||
output = OutputParameters(dirname=dirname+"/surface_fluxes", | ||
output = OutputParameters(dirname=dirname+"/wind_drag", | ||
There was a problem hiding this comment. Choose a reason for hiding this commentThe reason will be displayed to describe this comment to others. Learn more. Thanks for fixing this! |
||
dumpfreq=1, | ||
dumplist=['u']) | ||
io = IO(domain, output) | ||
|
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
Since the wind is zero, do we need a function prescribing the velocity? Or could we just set the initial velocity to zero?
There was a problem hiding this comment.
Choose a reason for hiding this comment
The reason will be displayed to describe this comment to others. Learn more.
I can't find a way to get the prescribed_transporting_velocity to work without giving it a function.