-
Notifications
You must be signed in to change notification settings - Fork 11
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
Allow reference profiles to be updated #542
Open
tommbendall
wants to merge
8
commits into
main
Choose a base branch
from
TBendall/update_ref_profile
base: main
Could not load branches
Branch not found: {{ refName }}
Loading
Could not load tags
Nothing to show
Loading
Are you sure you want to change the base?
Some commits from the old base branch may be removed from the timeline,
and old review comments may become outdated.
Open
Changes from 4 commits
Commits
Show all changes
8 commits
Select commit
Hold shift + click to select a range
f516f68
allow reference profiles to be updated
tommbendall cd133ba
get updating reference profile working with the help of the checkpoin…
tommbendall 9ab9337
Merge branch 'main' into TBendall/update_ref_profile
tommbendall 2a4dd5d
Merge branch 'main' into TBendall/update_ref_profile
tommbendall f4b43fe
Merge branch 'main' into TBendall/update_ref_profile
jshipton 6418788
used named tuple to encapsulate timing data
tommbendall 028fae6
Merge branch 'TBendall/update_ref_profile' of github.com:firedrakepro…
tommbendall a71d25c
Merge branch 'main' into TBendall/update_ref_profile
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 |
---|---|---|
|
@@ -35,7 +35,8 @@ def __init__(self, equation_set, io, transport_schemes, spatial_methods, | |
diffusion_schemes=None, physics_schemes=None, | ||
slow_physics_schemes=None, fast_physics_schemes=None, | ||
alpha=Constant(0.5), off_centred_u=False, | ||
num_outer=2, num_inner=2, accelerator=False): | ||
num_outer=2, num_inner=2, accelerator=False, | ||
reference_update_freq=None): | ||
|
||
""" | ||
Args: | ||
|
@@ -84,13 +85,24 @@ def __init__(self, equation_set, io, transport_schemes, spatial_methods, | |
implicit forcing (pressure gradient and Coriolis) terms, and the | ||
linear solve. Defaults to 2. Note that default used by the Met | ||
Office's ENDGame and GungHo models is 2. | ||
accelerator (bool, optional): Whether to zero non-wind implicit forcings | ||
for transport terms in order to speed up solver convergence | ||
accelerator (bool, optional): Whether to zero non-wind implicit | ||
forcings for transport terms in order to speed up solver | ||
convergence. Defaults to False. | ||
reference_update_freq (float, optional): frequency with which to | ||
update the reference profile with the n-th time level state | ||
fields. This variable corresponds to time in seconds, and | ||
setting this to zero will update the reference profiles every | ||
time step. Setting it to None turns off the update, and | ||
reference profiles will remain at their initial values. | ||
Defaults to None. | ||
""" | ||
|
||
self.num_outer = num_outer | ||
self.num_inner = num_inner | ||
self.alpha = alpha | ||
self.accelerator = accelerator | ||
self.reference_update_freq = reference_update_freq | ||
self.to_update_ref_profile = False | ||
|
||
# default is to not offcentre transporting velocity but if it | ||
# is offcentred then use the same value as alpha | ||
|
@@ -188,7 +200,6 @@ def __init__(self, equation_set, io, transport_schemes, spatial_methods, | |
self.linear_solver = linear_solver | ||
self.forcing = Forcing(equation_set, self.alpha) | ||
self.bcs = equation_set.bcs | ||
self.accelerator = accelerator | ||
|
||
def _apply_bcs(self): | ||
""" | ||
|
@@ -252,6 +263,24 @@ def copy_active_tracers(self, x_in, x_out): | |
for name in self.tracers_to_copy: | ||
x_out(name).assign(x_in(name)) | ||
|
||
def update_reference_profiles(self): | ||
""" | ||
Updates the reference profiles and if required also updates them in the | ||
linear solver. | ||
""" | ||
|
||
if self.reference_update_freq is not None: | ||
if float(self.t) + self.reference_update_freq > self.last_ref_update_time: | ||
self.equation.X_ref.assign(self.x.n(self.field_name)) | ||
self.last_ref_update_time = float(self.t) | ||
if hasattr(self.linear_solver, 'update_reference_profiles'): | ||
self.linear_solver.update_reference_profiles() | ||
|
||
elif self.to_update_ref_profile: | ||
if hasattr(self.linear_solver, 'update_reference_profiles'): | ||
self.linear_solver.update_reference_profiles() | ||
self.to_update_ref_profile = False | ||
|
||
def timestep(self): | ||
"""Defines the timestep""" | ||
xn = self.x.n | ||
|
@@ -264,13 +293,18 @@ def timestep(self): | |
xrhs_phys = self.xrhs_phys | ||
dy = self.dy | ||
|
||
# Update reference profiles -------------------------------------------- | ||
self.update_reference_profiles() | ||
|
||
# Slow physics --------------------------------------------------------- | ||
x_after_slow(self.field_name).assign(xn(self.field_name)) | ||
if len(self.slow_physics_schemes) > 0: | ||
with timed_stage("Slow physics"): | ||
logger.info('Semi-implicit Quasi Newton: Slow physics') | ||
for _, scheme in self.slow_physics_schemes: | ||
scheme.apply(x_after_slow(scheme.field_name), x_after_slow(scheme.field_name)) | ||
|
||
# Explict forcing ------------------------------------------------------ | ||
with timed_stage("Apply forcing terms"): | ||
logger.info('Semi-implicit Quasi Newton: Explicit forcing') | ||
# Put explicit forcing into xstar | ||
|
@@ -280,8 +314,10 @@ def timestep(self): | |
# the correct values | ||
xp(self.field_name).assign(xstar(self.field_name)) | ||
|
||
# OUTER ---------------------------------------------------------------- | ||
for outer in range(self.num_outer): | ||
|
||
# Transport -------------------------------------------------------- | ||
with timed_stage("Transport"): | ||
self.io.log_courant(self.fields, 'transporting_velocity', | ||
message=f'transporting velocity, outer iteration {outer}') | ||
|
@@ -290,6 +326,7 @@ def timestep(self): | |
# transports a field from xstar and puts result in xp | ||
scheme.apply(xp(name), xstar(name)) | ||
|
||
# Fast physics ----------------------------------------------------- | ||
x_after_fast(self.field_name).assign(xp(self.field_name)) | ||
if len(self.fast_physics_schemes) > 0: | ||
with timed_stage("Fast physics"): | ||
|
@@ -302,8 +339,7 @@ def timestep(self): | |
|
||
for inner in range(self.num_inner): | ||
|
||
# TODO: this is where to update the reference state | ||
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. Was this comment just wrong? |
||
|
||
# Implicit forcing --------------------------------------------- | ||
with timed_stage("Apply forcing terms"): | ||
logger.info(f'Semi-implicit Quasi Newton: Implicit forcing {(outer, inner)}') | ||
self.forcing.apply(xp, xnp1, xrhs, "implicit") | ||
|
@@ -314,6 +350,7 @@ def timestep(self): | |
xrhs -= xnp1(self.field_name) | ||
xrhs += xrhs_phys | ||
|
||
# Linear solve ------------------------------------------------- | ||
with timed_stage("Implicit solve"): | ||
logger.info(f'Semi-implicit Quasi Newton: Mixed solve {(outer, inner)}') | ||
self.linear_solver.solve(xrhs, dy) # solves linear system and places result in dy | ||
|
@@ -353,10 +390,18 @@ def run(self, t, tmax, pick_up=False): | |
pick_up: (bool): specify whether to pick_up from a previous run | ||
""" | ||
|
||
if not pick_up: | ||
if not pick_up and self.reference_update_freq is None: | ||
assert self.reference_profiles_initialised, \ | ||
'Reference profiles for must be initialised to use Semi-Implicit Timestepper' | ||
|
||
if not pick_up and self.reference_update_freq is not None: | ||
# Force reference profiles to be updated on first time step | ||
self.last_ref_update_time = float(t) - float(self.dt) | ||
|
||
elif not pick_up or (pick_up and self.reference_update_freq is None): | ||
# Indicate that linear solver profile needs updating | ||
self.to_update_ref_profile = True | ||
|
||
super().run(t, tmax, pick_up=pick_up) | ||
|
||
|
||
|
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
Oops, something went wrong.
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.
I feel slightly concerned about having a tuple of things that I might need to remember the order of! Could we use a named tuple?
Something like (completely untested!):
Then you can later access these like:
t = time_data.t
.