Using PyClaw’s solvers: Classic and SharpClaw¶
At present, PyClaw includes two types of solvers:
Classic: the original Clawpack algorithms, in 1/2/3D
SharpClaw: higher-order wave propagation using WENO reconstruction and Runge-Kutta integration, in 1/2D
Solver initialization takes one argument: a Riemann solver, usually from the Riemann repository. Typically, all that is needed to select a different solver is to specify it in the problem script, e.g.
>>> from clawpack import pyclaw
>>> from clawpack import riemann
>>> solver = pyclaw.ClawSolver2D(riemann.acoustics_2D)
for the 2D acoustics equations and the Classic Clawpack solver or
>>> solver = pyclaw.SharpClawSolver2D(riemann.acoustics_2D)
for the SharpClaw solver. Most of the applications distributed with PyClaw are set up to use either solver, depending on the value of the command line option solver_type, which should be set to classic or sharpclaw.
Typically, for a given grid resolution, the SharpClaw solvers are more accurate but also more computationally expensive. For typical problems involving shocks, the Classic solvers are recommended. For problems involving high-frequency waves, turbulence, or smooth solutions, the SharpClaw solvers may give more accurate solutions at less cost. This is an active area of research and you may wish to experiment with both solvers.
Future plans include incorporation of finite-difference and discontinuous Galerkin solvers.
Key differences between the Classic and SharpClaw solvers are:
The source term routine for the Classic solver should return the integral of the source term over a step, while the source term routine for SharpClaw should return the instantaneous value of the source term. For Classic, the source term function is set using solver.step_source, while for SharpClaw it is set using solver.dq_src. The shock-bubble interaction example shows how to use each of these.
The solvers have different options. For a list of options and possible values, see the documentation of the
ClawSolver
andSharpClawSolver
classes.
Pyclaw Classic Clawpack Solvers¶
The pyclaw classic clawpack solvers are a collection of solvers that represent the functionality of the older versions of clawpack. It comes in two forms, a pure python version and a python wrapping of the fortran libraries. All of the solvers available provide the same basic interface and provide the same options as the old versions of clawpack. The superclass solvers are not meant to be used separately but there to provide common routines for all the Clawpack solvers. Please refer to each of the inherited classes for more info about the methods and attributes they provide each class. .. The inheritance structure is:
- Example:
This is a simple example of how to instantiate and evolve a solution to a later time \(\text{t_end}\) using the linearized 1d acoustics Riemann solver
>>> from clawpack import pyclaw
>>> solver = pyclaw.ClawSolver1D() # Instantiate a default, 1d solver
>>> solver.limiters = pyclaw.limiters.tvd.vanleer # Use the van Leer limiter
>>> solver.dt = 0.0001 # Set the initial time step
>>> solver.max_steps = 500 # Set the maximum number of time steps
>>> solver.evolve_to_time(solution,t_end) # Evolve the solution to t_end
pyclaw.classic.solver
¶
- class clawpack.pyclaw.classic.solver.ClawSolver(riemann_solver=None, claw_package=None)¶
Generic classic Clawpack solver
All Clawpack solvers inherit from this base class.
- mthlim¶
Limiter(s) to be used. Specified either as one value or a list. If one value, the specified limiter is used for all wave families. If a list, the specified values indicate which limiter to apply to each wave family. Take a look at pyclaw.limiters.tvd for an enumeration.
Default = limiters.tvd.minmod
- order¶
Order of the solver, either 1 for first order (i.e., Godunov’s method) or 2 for second order (Lax-Wendroff-LeVeque).
Default = 2
- source_split¶
Which source splitting method to use: 1 for first order Godunov splitting and 2 for second order Strang splitting.
Default = 1
- fwave¶
Whether to split the flux jump (rather than the jump in Q) into waves; requires that the Riemann solver performs the splitting.
Default = False
- step_source¶
Handle for function that evaluates the source term. The required signature for this function is:
def step_source(solver,state,dt)
- kernel_language¶
Specifies whether to use wrapped Fortran routines (‘Fortran’) or pure Python (‘Python’).
Default = 'Fortran'
.
- verbosity¶
The level of detail of logged messages from the Fortran solver.
Default = 0
.
- setup(solution)¶
Perform essential solver setup. This routine must be called before solver.step() may be called.
- step(solution, take_one_step, tstart, tend)¶
Evolve solution one time step
The elements of the algorithm for taking one step are:
Pick a step size as specified by the base solver attribute
get_dt()
A half step on the source term
step_source()
if Strang splitting is being used (source_split
= 2)A step on the homogeneous problem \(q_t + f(q)_x = 0\) is taken
A second half step or a full step is taken on the source term
step_source()
depending on whether Strang splitting was used (source_split
= 2) or Godunov splitting (source_split
= 1)
This routine is called from the method evolve_to_time defined in the pyclaw.solver.Solver superclass.
- Input:
solution - (
Solution
) solution to be evolved
- Output:
(bool) - True if full step succeeded, False otherwise
- step_hyperbolic(solution)¶
Take one homogeneous step on the solution.
This is a dummy routine and must be overridden.
Change to Custom BC Function Signatures¶
To allow better access to aux array data in the boundary condition functions both the qbc and auxbc arrays are now passed to the custom boundary condition functions. The new signature is
- def my_custom_BC(state, dim, t, qbc, auxbc, num_ghost):
…
and should be adopted as soon as possible. The old signature
- def my_custom_BC(state, dim, t, bc_array, num_ghost):
…
can still be used but a warning will be issued and the old signature will not be supported when version 6.0 is released. This addition is available in versions > 5.2.0.