37 New Zealand#

Authors: Kidus Teshome, Cian Wilson

Steady state 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, 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_problems.sz_params import allsz_params
from sz_problems.sz_slab import create_slab, plot_slab
from sz_problems.sz_geometry import create_sz_geometry
from sz_problems.sz_steady_dislcreep import SteadyDislSubductionProblem
import numpy as np
import dolfinx as df
import pyvista as pv
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 of the model.

Resolution

By default the resolution is low to allow for a quick runtime and smaller website size. If sufficient computational resources are available set a lower resscale to get higher resolutions and results with sufficient accuracy.

name = "37_New_Zealand"
resscale = 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))
37_New_Zealand:
Key                  Value     
-------------------------------------------------------------------------------------
coast_distance       185
extra_width          18.5
lc_depth             30
io_depth             168
uc_depth             15
dirname              37_New_Zealand
As                   40.0
qs                   0.065
A                    100
sztype               continental
Vs                   30.4
xs                   [0, 41.5, 89.5, 167.0, 180.2, 183.2, 200.0, 259.0, 281.5]
ys                   [-6, -15.0, -30.0, -70.0, -80.0, -82.5, -100.0, -200.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 SteadyDislSubductionProblem class.

SteadyDislSubductionProblem?

Setup#

Setup a slab.

slab = create_slab(szdict['xs'], szdict['ys'], resscale, szdict['lc_depth'])
_ = plot_slab(slab)
../../../_images/bd353c37bde29700cc2f029b9cfd67d1e5600e5836041fcd19e10b98b15929c1.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/0e936a0676a5b49d9ec7f6901a1a48b5865654441c09fe0f913935c3e7f8ddda.png

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

sz = SteadyDislSubductionProblem(geom, **szdict)

Solve#

Solve using a dislocation creep rheology and assuming a steady state.

sz.solve()
Iteration   Residual     Relative Residual
------------------------------------------
0           6143.77      1           
1           4082.68      0.664524    
2           174.957      0.0284772   
3           75.3323      0.0122616   
4           30.0221      0.0048866   
5           14.9905      0.00243996  
6           9.64427      0.00156976  
7           7.18037      0.00116872  
8           5.59524      0.000910718 
9           4.41493      0.000718604 
10          3.49684      0.000569169 
11          2.77684      0.000451976 
12          2.20977      0.000359677 
13          1.75747      0.000286057 
14          1.39659      0.000227319 
15          1.10945      0.000180581 
16          0.880359     0.000143293 
17          0.697555     0.000113539 
18          0.551861     8.98245e-05 
19          0.43596      7.09597e-05 
20          0.343951     5.59838e-05 
21          0.271056     4.41188e-05 
22          0.213406     3.47354e-05 
23          0.167887     2.73265e-05 
24          0.131997     2.14847e-05 
25          0.10373      1.68837e-05 
26          0.0814865    1.32633e-05 
27          0.0639985    1.04168e-05 
28          0.0502577    8.18028e-06 
29          0.0394662    6.42378e-06 
30          0.0309934    5.04468e-06 
31          0.0243421    3.96209e-06 

Plot#

Plot the solution.

plotter = pv.Plotter()
utils.plot.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.plot_vector_glyphs(sz.vw_i, plotter=plotter, gather=True, factor=0.1, color='k', scale=utils.mps_to_mmpyr(sz.v0))
utils.plot.plot_vector_glyphs(sz.vs_i, plotter=plotter, gather=True, factor=0.1, color='k', scale=utils.mps_to_mmpyr(sz.v0))
utils.plot.plot_geometry(geom, plotter=plotter, color='green', width=2)
utils.plot.plot_couplingdepth(slab, plotter=plotter, render_points_as_spheres=True, point_size=10.0, color='green')
utils.plot.plot_show(plotter)
utils.plot.plot_save(plotter, output_folder / "{}_ss_solution_resscale_{:.2f}.png".format(name, resscale))

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

filename = output_folder / "{}_ss_solution_resscale_{:.2f}.bp".format(name, resscale)
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
zipfilename = filename.with_suffix(".zip")
!zip -r $zipfilename $filename
  adding: output/37_New_Zealand_ss_solution_resscale_3.00.bp/ (stored 0%)
  adding: output/37_New_Zealand_ss_solution_resscale_3.00.bp/mmd.0 (deflated 64%)
  adding: output/37_New_Zealand_ss_solution_resscale_3.00.bp/profiling.json (deflated 62%)
  adding: output/37_New_Zealand_ss_solution_resscale_3.00.bp/md.idx (deflated 59%)
  adding: output/37_New_Zealand_ss_solution_resscale_3.00.bp/md.0 (deflated 75%)
  adding: output/37_New_Zealand_ss_solution_resscale_3.00.bp/data.0 (deflated 63%)

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, 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_ss', szdict['dirname']+'_minres_2.00.vtu')
tffilepath = os.path.join(basedir, os.path.pardir, os.path.pardir, os.path.pardir, 'data')
with zipfile.ZipFile(zipfilename, 'r') as z:
    z.extract(tffilename, path=tffilepath)
fxgrid = utils.plot.grids_scalar(sz.T_i)[0]

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

diffgrid = utils.plot.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.plot_geometry(geom, plotter=plotter_diff, color='green', width=2)
utils.plot.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.3415587338032692