Step 5: Static Coseismic Slip#

Features

  • Triangular cells

  • pylith.meshio.MeshIOPetsc

  • pylith.problems.TimeDependent

  • pylith.bc.DirichletTimeDependent

  • spatialdata.spatialdb.SimpleDB

  • spatialdata.spatialdb.ZeroDB

  • pylith.meshio.OutputSolnBoundary

  • pylith.meshio.DataWriterHDF5

  • Static simulation

  • pylith.materials.Elasticity

  • pylith.materials.IsotropicLinearElasticity

  • pylith.faults.FaultCohesiveKin

  • pylith.faults.KinSrcStep

  • spatialdata.spatialdb.UniformDB

Simulation parameters#

This example involves a static simulation that solves for the deformation from prescribed coseismic slip on the main fault. We specify 2 meters of reverse slip. Fig. 79 shows the boundary conditions on the domain. The parameters specific to this example are in step05_onefault.cfg.

Fig. 79 Boundary conditions for static coseismic slip on the main fault. We prescribe 2 meters of reverse slip with roller boundary conditions on the lateral sides and bottom of the domain.#

Important

In 2D simulations slip is specified in terms of opening and left-lateral components. This provides a consistent, unique sense of slip that is independent of the fault orientation. For our geometry in this example, right lateral slip corresponds to reverse slip on the dipping fault.

Important

The main fault contains one end that is buried within the domain. When PyLith inserts cohesive cells into a mesh with buried edges (in this case a point), we must identify these buried edges so that PyLith properly adjusts the topology along these edges.

We adjust the solution field to include both displacement and the Lagrange multiplier associated with the fault. For uniform prescribed slip we use a UniformDB.

Listing 114 Parameters for earthquake rupture on the main reverse fault in Step 5.#
[pylithapp.problem]
solution = pylith.problems.SolnDispLagrange

[pylithapp.problem.interfaces.fault]
label = fault
label_value = 20
edge = fault_end
edge_value = 21
observers.observer.data_fields = [slip]

[pylithapp.problem.interfaces.fault.eq_ruptures.rupture]
db_auxiliary_field = spatialdata.spatialdb.UniformDB
db_auxiliary_field.description = Fault rupture auxiliary field spatial database
db_auxiliary_field.values = [initiation_time, final_slip_left_lateral, final_slip_opening]
db_auxiliary_field.data = [0.0*s, -2.0*m, 0.0*m]

Running the simulation#

Listing 115 Run Step 5 simulation#
$ pylith step05_onefault.cfg

# The output should look something like the following.
 >> /software/unix/py39-venv/pylith-debug/lib/python3.9/site-packages/pylith/meshio/MeshIOObj.py:44:read
 -- meshiopetsc(info)
 -- Reading finite-element mesh
 >> /src/cig/pylith/libsrc/pylith/meshio/MeshIO.cc:94:void pylith::meshio::MeshIO::read(pylith::topology::Mesh*)
 -- meshiopetsc(info)
 -- Component 'reader': Domain bounding box:
    (-100000, 100000)
    (-100000, 0)

# -- many lines omitted --

 >> /software/unix/py39-venv/pylith-debug/lib/python3.9/site-packages/pylith/problems/TimeDependent.py:139:run
 -- timedependent(info)
 -- Solving problem.
0 TS dt 0.01 time 0.
    0 SNES Function norm 2.129295960330e-02
    Linear solve converged due to CONVERGED_ATOL iterations 43
    1 SNES Function norm 1.577267844919e-12
  Nonlinear solve converged due to CONVERGED_FNORM_ABS iterations 1
1 TS dt 0.01 time 0.01
 >> /software/unix/py39-venv/pylith-debug/lib/python3.9/site-packages/pylith/problems/Problem.py:201:finalize
 -- timedependent(info)
 -- Finalizing problem.

From the end of the output written to the terminal window, we see that the linear solver converged in 43 iterations and met the absolute convergence tolerance (ksp_atol). As we expect for this linear problem, the nonlinear solver converged in 1 iteration.

Visualizing the results#

The output directory contains the simulation output. Each “observer” writes its own set of files, so the solution over the domain is in one set of files, the boundary condition information is in another set of files, and the material information is in yet another set of files. The HDF5 (.h5) files contain the mesh geometry and topology information along with the solution fields. The Xdmf (.xmf) files contain metadata that allow visualization tools like ParaView to know where to find the information in the HDF5 files. To visualize the data using ParaView or Visit, load the Xdmf files.

In Fig. 80 we use ParaView to visualize the x displacement field using the viz/plot_dispwarp.py Python script. First, we start ParaView from the examples/reverse-2d directory. Before running the viz/plot_dispwarp.py Python script as described in ParaView Python Scripts, we set the simulation name in the ParaView Python Shell.

Listing 116 Set the simulation and exaggeration in the ParaView Python Shell.#
>>> SIM = "step05_onefault"
>>> FIELD_COMPONENT = "X"
Solution for Step 5. The colors indicate the magnitude of the x displacement, and the deformation is exaggerated by a factor of 1000.

Fig. 80 Solution for Step 5. The colors of the shaded surface indicate the magnitude of the x displacement, and the deformation is exaggerated by a factor of 1000. The undeformed configuration is show by the gray wireframe.#