06 GuatElSal#

Authors: Cameron Seebeck, Cian Wilson

Time-dependent implementation#

Preamble#

Set some path information.

import sys, os
basedir = ''
if "__file__" in globals(): basedir = os.path.dirname(__file__)
sys.path.append(os.path.join(basedir, os.path.pardir))
sys.path.append(os.path.join(basedir, os.path.pardir, os.path.pardir, 'python'))

Loading everything we need from sz_problem and also set our default plotting and output preferences.

import utils
from sz_base import allsz_params
from sz_slab import create_slab, plot_slab
from sz_geometry import create_sz_geometry
from sz_problem import SubductionProblem
import numpy as np
import dolfinx as df
import pyvista as pv
if __name__ == "__main__" and "__file__" in globals():
    pv.OFF_SCREEN = True
import pathlib
output_folder = pathlib.Path(os.path.join(basedir, "output"))
output_folder.mkdir(exist_ok=True, parents=True)
import hashlib
import zipfile
import requests

Parameters#

We first select the name and resolution scale, resscale and target Courant number cfl of the model.

Resolution

By default the resolution (both spatial and temporal) is low to allow for a quick runtime and smaller website size. If sufficient computational resources are available set a lower resscale and a lower cfl to get higher spatial and temporal resolutions respectively. This is necessary to get results with sufficient accuracy for scientific interpretation.

name = "06_GuatElSal"
resscale = 3.0
cfl      = 3.0

Then load the remaining parameters from the global suite.

szdict = allsz_params[name]
print("{}:".format(name))
print("{:<20} {:<10}".format('Key','Value'))
print("-"*85)
for k, v in allsz_params[name].items():
    if v is not None and k not in ['z0', 'z15']: print("{:<20} {}".format(k, v))
06_GuatElSal:
Key                  Value     
-------------------------------------------------------------------------------------
coast_distance       155
extra_width          20
lc_depth             45
io_depth             177
uc_depth             15
dirname              06_GuatElSal
As                   40.0
qs                   0.065
A                    17.4
sztype               continental
Vs                   67.0
xs                   [0, 41.6, 104.0, 129.0, 152.5, 154.9, 169.0, 204.0, 245.0]
ys                   [-6, -15.0, -45.0, -60.0, -80.0, -82.5, -100.0, -160.0, -240.0]

Any of these can be modified in the dictionary.

Several additional parameters can be modified, for details see the documentation for the SubductionProblem class.

if __name__ == "__main__" and "__file__" not in globals():
    SubductionProblem?

The if __name__ == "__main__" and "__file__" not in globals(): logic above is only necessary to make sure that this only runs in the Jupyter notebook version of this code and not the python version. It is not generally necessary when getting the docstring of a function or class in Jupyter.

Setup#

Setup a slab.

slab = create_slab(szdict['xs'], szdict['ys'], resscale, szdict['lc_depth'])
_ = plot_slab(slab)
../../_images/2d3c0260cecc1dbf84a3b8c29ffa08e7b45cba8571c8f28dda091fd1cf57349f.png

Create the subduction zome geometry around the slab.

geom = create_sz_geometry(slab, resscale, szdict['sztype'], szdict['io_depth'], szdict['extra_width'], 
                             szdict['coast_distance'], szdict['lc_depth'], szdict['uc_depth'])
_ = geom.plot()
../../_images/f673affa3ccee4989c15dab1aa074003dc8eb2821368ce7daa59ac026866e64c.png

Finally, declare the SubductionZone problem class using the dictionary of parameters.

sz = SubductionProblem(geom, **szdict)

Solve#

Let’s first set up a gif of the temperature.

fps = 5
plotter = pv.Plotter(notebook=False, off_screen=True)
utils.plot_scalar(sz.T_i, plotter=plotter, scale=sz.T0, gather=True, cmap='coolwarm', clim=[0.0, sz.Tm*sz.T0], scalar_bar_args={'title': 'Temperature (deg C)', 'bold':True})
utils.plot_geometry(geom, plotter=plotter, color='green', width=2)
utils.plot_couplingdepth(slab, plotter=plotter, render_points_as_spheres=True, point_size=10.0, color='green')
plotter.open_gif( str(output_folder / "{}_td_solution_resscale_{:.2f}_cfl_{:.2f}.gif".format(name, resscale, cfl,)), fps=fps)

Solve using a dislocation creep rheology.

# Select the timestep based on the approximate target Courant number
dt = cfl*resscale/szdict['Vs']
# Reduce the timestep to get an integer number of timesteps
dt = szdict['As']/np.ceil(szdict['As']/dt)
sz.solve_timedependent_dislocationcreep(szdict['As'], dt, theta=0.5, rtol=1.e-1, verbosity=1, plotter=plotter)
plotter.close()
Entering timeloop with 298 steps (dt = 0.134228 Myr, final time = 40 Myr)
Finished timeloop after 298 steps (final time = 40 Myr)

Plot#

Plot the solution.

plotter = pv.Plotter()
utils.plot_scalar(sz.T_i, plotter=plotter, scale=sz.T0, gather=True, cmap='coolwarm', scalar_bar_args={'title': 'Temperature (deg C)', 'bold':True})
utils.plot_vector_glyphs(sz.vw_i, plotter=plotter, gather=True, factor=0.1, color='k', scale=utils.mps_to_mmpyr(sz.v0))
utils.plot_vector_glyphs(sz.vs_i, plotter=plotter, gather=True, factor=0.1, color='k', scale=utils.mps_to_mmpyr(sz.v0))
utils.plot_geometry(geom, plotter=plotter, color='green', width=2)
utils.plot_couplingdepth(slab, plotter=plotter, render_points_as_spheres=True, point_size=10.0, color='green')
utils.plot_show(plotter)
utils.plot_save(plotter, output_folder / "{}_td_solution_resscale_{:.2f}_cfl_{:.2f}.png".format(name, resscale, cfl,))

Save it to disk so that it can be examined with other visualization software (e.g. Paraview).

filename = output_folder / "{}_td_solution_resscale_{:.2f}_cfl_{:.2f}.bp".format(name, resscale, cfl,)
with df.io.VTXWriter(sz.mesh.comm, filename, [sz.T_i, sz.vs_i, sz.vw_i]) as vtx:
    vtx.write(0.0)
# zip the .bp folder so that it can be downloaded from Jupyter lab
if __name__ == "__main__" and "__file__" not in globals():
    zipfilename = filename.with_suffix(".zip")
    !zip -r $zipfilename $filename
  adding: output/06_GuatElSal_td_solution_resscale_3.00_cfl_3.00.bp/ (stored 0%)
  adding: output/06_GuatElSal_td_solution_resscale_3.00_cfl_3.00.bp/md.idx (deflated 59%)
  adding: output/06_GuatElSal_td_solution_resscale_3.00_cfl_3.00.bp/profiling.json (deflated 61%)
  adding: output/06_GuatElSal_td_solution_resscale_3.00_cfl_3.00.bp/mmd.0 (deflated 64%)
  adding: output/06_GuatElSal_td_solution_resscale_3.00_cfl_3.00.bp/data.0 (deflated 64%)
  adding: output/06_GuatElSal_td_solution_resscale_3.00_cfl_3.00.bp/md.0 (deflated 76%)

Comparison#

Compare to the published result from Wilson & van Keken, PEPS, 2023 (II) and van Keken & Wilson, PEPS, 2023 (III). The original models used in these papers are also available as open-source repositories on github and zenodo.

First download the minimal necessary data from zenodo and check it is the right version.

zipfilename = pathlib.Path(os.path.join(basedir, os.path.pardir, os.path.pardir, "data", "vankeken_wilson_peps_2023_TF_lowres_minimal.zip"))
if not zipfilename.is_file():
    zipfileurl = 'https://zenodo.org/records/13234021/files/vankeken_wilson_peps_2023_TF_lowres_minimal.zip'
    r = requests.get(zipfileurl, allow_redirects=True)
    open(zipfilename, 'wb').write(r.content)
assert hashlib.md5(open(zipfilename, 'rb').read()).hexdigest() == 'a8eca6220f9bee091e41a680d502fe0d'
tffilename = os.path.join('vankeken_wilson_peps_2023_TF_lowres_minimal', 'sz_suite_td', szdict['dirname']+'_minres_2.00_cfl_2.00.vtu')
tffilepath = os.path.join(os.pardir, os.pardir, 'data')
with zipfile.ZipFile(zipfilename, 'r') as z:
    z.extract(tffilename, path=tffilepath)
fxgrid = utils.grids_scalar(sz.T_i)[0]

tfgrid = pv.get_reader(os.path.join(tffilepath, tffilename)).read()

diffgrid = utils.pv_diff(fxgrid, tfgrid, field_name_map={'T':'Temperature::PotentialTemperature'}, pass_point_data=True)
diffgrid.set_active_scalars('T')
plotter_diff = pv.Plotter()
clim = None
plotter_diff.add_mesh(diffgrid, cmap='coolwarm', clim=clim, scalar_bar_args={'title': 'Temperature Difference (deg C)', 'bold':True})
utils.plot_geometry(geom, plotter=plotter_diff, color='green', width=2)
utils.plot_couplingdepth(slab, plotter=plotter_diff, render_points_as_spheres=True, point_size=5.0, color='green')
plotter_diff.enable_parallel_projection()
plotter_diff.view_xy()
plotter_diff.show()
integrated_data = diffgrid.integrate_data()
error = integrated_data['T'][0]/integrated_data['Area'][0]
print("Average error = {}".format(error,))
assert np.abs(error) < 5
Average error = 0.47323597510459886

Finish up#

Convert this notebook to a python script (making sure to save first)

if __name__ == "__main__" and "__file__" not in globals():
    from ipylab import JupyterFrontEnd
    app = JupyterFrontEnd()
    app.commands.execute('docmanager:save')
    !jupyter nbconvert --NbConvertApp.export_format=script --ClearOutputPreprocessor.enabled=True 06_GuatElSal.ipynb
[NbConvertApp] Converting notebook 06_GuatElSal.ipynb to script
[NbConvertApp] Writing 8120 bytes to 06_GuatElSal.py