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 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 t_end using the linearized 1d acoustics Riemann solver:

>>> from pyclaw.evolve.clawpack import ClawSolver1D
>>> solver = ClawSolver1D()    # Instantiate a default, 1d solver
>>> solver.mthlim = [3,3]      # 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

Note

The fortran wrapped classes are not included in this release.

pyclaw.evolve.clawpack

Module containg the classic Clawpack solvers

This module contains the pure and wrapped classic clawpack solvers. All clawpack solvers inherit from the ClawSolver superclass which in turn inherits from the Solver superclass. As such, the only solver classes that should be directly used should be the dimensionally dependent ones such as ClawSolver1D.

Authors:Kyle T. Mandli (2008-09-11) Initial version
class pyclaw.evolve.clawpack.ClawSolver(data=None)

Bases: pyclaw.evolve.solver.Solver

Generic classic Clawpack solver

All Clawpack solvers inherit from this base class.

mthlim

Limiter to be used on each wave. Default = [1]

order

Order of the solver, either 1 for first order or 2 for second order corrections. Default = 2

src_split

Whether to use a source splitting method, 0 for none, 1 for first order Godunov splitting and 2 for second order Strang splitting. Default = 0

fwave

Whether to split the flux into waves, requires that the Riemann solver performs the splitting. Default = False

src

Source term function. Default is the stub function.

start_step

Function called before each time step is taken. Default is the stub function

Initialization :
 
Input:
  • data - (Data) Data object, the solver will look for the named variables to instantiate itself.
Output:
Version :1.0 (2009-06-01)
homogeneous_step(solutions)

Take one homogeneous step on the solutions

This is a dummy routine and must be overridden.

list_riemann_solvers()

List available Riemann solvers

This routine returns a list of available Riemann solvers which is constructed in the Riemann solver package (Riemann Solver Package). In this case it lists all Riemann solvers.

Output :

Note

These Riemann solvers are currently only accessible to the python time stepping routines.

set_riemann_solver(solver_name)

Assigns the library solver solver_name as the Riemann solver.

Input :
  • solver_name - (string) Name of the solver to be used, raises a NameError if the solver does not exist.
setup()

Called before any set of time steps.

This routine will be called once before the solver is used via the Controller. In the case of ClawSolver we make sure that the mthlim is a list.

step(solutions)

Evolve solutions one time step

This routine encodes the generic order in a full time step in this order:

  1. The start_step() function is called
  2. A half step on the source term src() if Strang splitting is being used (src_split = 2)
  3. A step on the homogeneous problem q_t + f(q)_x = 0 is taken
  4. A second half step or a full step is taken on the source term src() depending on whether Strang splitting was used (src_split = 2) or Godunov splitting (src_split = 1)

This routine is called from the method evolve_to_time defined in the pyclaw.evolve.solver.Solver superclass.

Input :
  • solutions - (Solution) Dictionary of solutions to be evolved
Output :
  • (bool) - True if full step succeeded, False otherwise
class pyclaw.evolve.clawpack.ClawSolver1D(data=None)

Bases: pyclaw.evolve.clawpack.ClawSolver

Clawpack evolution routine in 1D

This class represents the 1d clawpack solver on a single grid. Note that there are routines here for interfacing with the fortran time stepping routines and the python time stepping routines. The ones used are dependent on the argument given to the initialization of the solver (defaults to python).

rp

Riemann solver function.

Initialization :
 
Input:
  • data - (Data) An instance of a Data object whose parameters can be used to initialize this solver
Output:
Authors :Kyle T. Mandli (2008-09-11) Initial version
homogeneous_step(solutions)

Take one time step on the homogeneous hyperbolic system

Takes one time step of size dt on the hyperbolic system defined in the appropriate Riemann solver rp.

Input :
  • solutions - (Solution) Solution that will be evolved
Version :

1.0 (2009-07-01)

list_riemann_solvers()

List available Riemann solvers

This routine returns a list of available Riemann solvers which is constructed in the Riemann solver package (_pyclaw_rp). In this case it lists only the 1D Riemann solvers.

Output :

Note

These Riemann solvers are currently only accessible to the python time stepping routines.

set_riemann_solver(solver_name)

Assigns the library solver solver_name as the Riemann solver.

Input :
  • solver_name - (string) Name of the solver to be used, raises a NameError if the solver does not exist.
pyclaw.evolve.clawpack.src(solver, solutions, t, dt)

Dummy routine called to calculate a source term

Replace this routine if you want to include a source term.

pyclaw.evolve.clawpack.start_step(solver, solutions)

Dummy routine called before each step

Replace this routine if you want to do something before each time step.

Table Of Contents

Previous topic

Pyclaw Evolve Package

Next topic

Pyclaw Limiters

This Page