Preconditioners

Preconditioners

This section describes the different catagories of preconditioners, including the API and the requirements for creating new preconditioners. Note that whenever a linear solver is used, the Linear Solver interface should be used rather than the Linear Operator interface described here.

Type hierarchy

When creating a new preconditioner, it is very important to make it inherit from the proper supertype.

Abstract supertype of all preconditioner types. Preconditioners can be used with iterative methods only. When using a direct method, PCNone should be used.

The purpose of this type is to provide a consistent interface for different types of preconditioners. In particular, it provides control over when the preconditioner is recomputed.

Users should implement new preconditioners using composition with one of the preconditioners defined here, namely PetscMatPC or PetscMatFreePC, ie. they should define a new subtype of AbstractPC that has either PetscMatPC or PetscMatFrePC as a field. This allows calling the existing functions for these types (which compute a preconditioner for the Jacobian of the physics) and modifying the preconditioner as needed. User defined preconditioners should subtype either AbstractPetscMatPC of AbstractPetscMatFreePC.

Fields

Note that arbitrarily deep nesting of preconditioners is allowed. The pc_inner field can be one of PCNone, PetscMatPC, or PetscMatFreePC for a non-nested preconditioner, or some other AbstractPC for a nested preconditioner.

source

Abstract supertype of all Petsc matrix-explicit preconditioners

source

Abstract supertype of all Petsc matrix-free preconditioners.

source

Alias for any kind of Petsc PC (matrix-explicit or matrix-free)

source

API

Every preconditioner supports the following functions. When defining a new preconditioner, some of the functions must be defined for the new AbstractPC type, while others are defined automatically based on the supertype of the preconditioner.

This function calculates the preconditioner. Every implementation of AbstractPC should extend this function with a new method

When creating new preconditioners, this function should generally be called first on pc_inner, and modifications to the Jacobian should be made subsequently.

Inputs

  • pc: the AbstractPC implementation

  • mesh

  • sbp

  • eqn

  • opts

  • ctx_residual: the ctx required by physicsRhs

  • t: current time

Implementation Notes: For matrix-explicit preconditioners, this might not actually calculate the preconditioner. Rather, it calculates the matrix the preconditioner is based on, and the solve function calcultes the preconditioner from it. Nevertheless, it supports the proper semantics for when the PC is updated even when the PC matrix and LinearOperator matrix are the same (as long as the solve function is called regularly).

source

Applies the preconditioner, ie. x = inv(Ap)*b, where Ap is the approximation to the matrix A. Note that Ap itself may not be available for some preconditioners, hence there is only an API for applying inv(Ap), not Ap itself.

Matrix-free preconditioners need to extend this function with a new method,
matrix-explicit preconditioners do not.

Inputs

  • pc: the AbstractPC implementation.

  • mesh

  • sbp

  • eqn: this argument should generally not be used because all the solution related data should be stored in the pc object by calcPC

  • opts

  • t: current time

  • b: a AbstractVector representing the local part of the solution (ie eqn.q_vec)

Inputs/Outputs

  • x: AbstractVector updated with results (same size as b) (do not overwrite)

source

Applies the transpose of the preconditioner, ie. x = inv(Ap).'*b. Similar to applyPC, see that function for details.

Note that not every preconditioning method supports this.

source

This function returns the underlying preconditioner object, ie. PCNone, PetscMatPC, or PetscMatFreePC.

Note that arbitrarily deep nesting of preconditioners is allowed. Users do not have to implement as long as the nested preconditioner is stored in a field called pc_inner.

For matrix-explicit preconditioners, this function is useful for getting the PetscMatPC object, which contains the preconditioning Jacobian matrix.

Inputs

source
Utils.freeMethod.

This function frees any memory belonging to external libraries. Users must call this function when they are finished with an AbstractPC object.

Users do not have to define this function for their AbstractPC types.

Inputs

  • pc: the AbstractPC object

source

Concrete PC Types

The following types are the "Base" PC types refered to by getBasePC. Every preconditioner must contain one of these (either directly, in the pc_inner field described by AbstractPC, of inside another preconditioner).

Preconditioner type for direct solve. Do not use with PetscLinearOperator.

Public Fields

  • none

source

Outer constructor

Inputs

  • mesh

  • sbp

  • eqn

  • opts

source

AbstractPC implementation for Petsc matrix-explicit preconditioners.

Public Fields

  • A: a PetscMat object used to calculate the preconditioner

source

Outer constructor

Inputs

  • mesh

  • sbp

  • eqn

  • opts

source

Type for managing a Petsc matrix-free preconditioner. Actual preconditioners can be built on top of this. Because matrix-free preconditioners do not have a particular form, this type requires the preconditioner to implement the entire AbstractPC interface. The benefit of using this type to build a preconditioner is that is automatically exposes the preconditioner to Petsc.

The calcPC function the user defines must call setPCCtx.

Public Fields

  • none

source

Outer constructor

Inputs

  • mesh

  • sbp

  • eqn

  • opts

source

PetscMatFreePC has an additional function not needed by other preconditioner types:

This function should be called by the users calcPC function. The arguments passed to this function are passed to applyPC during the next preconditioner application. See that function for a description of the arguments

Inputs

  • pc:

  • mesh

  • sbp

  • eqn

  • opts

  • ctx_residual

  • t

source

Implementing a New Preconditioner

Many of the API functions are automatically generated for matrix-explicit preconditioners. Matrix-free preconditioners need to define most of the PC interface themselves. A summary of the functions each preconditioner must implements is:

Matrix-explicit

Matrix-free

First Level PC

This example shows how to implement a new preconditioner that directly contains one of the Concrete PC Types described above.

type ExamplePetscMatPC <: AbstractPetscMatPC
  pc_inner::PetscMatPC  # the supertype of pc_inner and ExamplePetscMatPC
                        # must match
end

function ExamplePetscMatPC(mesh::AbstractMesh, sbp::AbstractSBP,
                    eqn::AbstractSolutionData, opts::Dict)

  pc_inner = PetscMatPC(mesh, sbp, eqn, opts)

  return ExamplePetscMatPC(pc_inner)
end


function calcPC(pc::ExamplePetscMatPC, mesh::AbstractMesh, sbp::AbstractSBP,
                eqn::AbstractSolutionData, opts::Dict, ctx_residual, t)


  calcPC(pc.pc_inner)

  pc2 = getBasePC(pc)  # because pc_inner is the base pc, this returns
                       # pc.pc_inner
  # call set_values1!(pc2.Ap, ...) to put values into the preconditioning matrix

  return nothing
end

Because calcPC is defined and the base PC is one of the concrete PC types, namely PetscMatPC, the functions applyPC and applyPC are automatically defined for ExamplePetscMatPC. The new preconditioner is now fully usable.

Multiple Levels of Composition

The same structure as shown in the previous structure can be used to build a new preconditioner out of any existing preconditioner, not only the Concrete PC Types.

type OuterPetscMatPC <: AbstractPetscMatPC
  pc_inner::ExamplePetscMatPC  # use the PC define in the previous example
end

function OuterPetscMatPC(mesh::AbstractMesh, sbp::AbstractSBP,
                         eqn::AbstractSolutionData, opts::Dict)

  pc_inner = ExamplePetscMatPC(mesh, sbp, eqn, opts)

  return OuterPetscMatPC(pc_inner)
end


function calcPC(pc::OuterPetscMatPC, mesh::AbstractMesh, sbp::AbstractSBP,
                eqn::AbstractSolutionData, opts::Dict, ctx_residual, t)


  calcPC(pc.pc_inner)  # always call the inner PC calc function first

  pc2 = getBasePC(pc)  # in this case, pc2 is the PetscMatPC inside 
                       # the ExamplePetscMatPC
  # call set_values1!(pc2.Ap, ...) to put values into the preconditioning matrix  # because calcPC(pc.pc_inner) was called already, the values here should
  # be added to those already in pc2.Ap

  return nothing
end

As before, this preconditioner is now fully usable.

Matrix Free

An outline of how a matrix free preconditioner should be implemented is:

type ExamplePetscMatFreePC <: AbstractPetscMatFreePC
  pc_inner::PetscMatFreePC
  # other fields...
end

function ExamplePetscMatFreePC(mesh::AbstractMesh, sbp::AbstractSBP,
                         eqn::AbstractSolutionData, opts::Dict)

  pc_inner = PetscMatFreePC(mesh, sbp, eqn, opts)

  return ExamplePetscMatFreePC(pc_inner)
end

function calcPC(pc::ExamplePetscMatFreePC, mesh::AbstractMesh, sbp::AbstractSBP,
                eqn::AbstractSolutionData, opts::Dict, ctx_residual, t)

  # do any setup operations, storing data into the "other fields..." of the pc

  # these arguments will be passed into applyPC() the next time it is called
  # by Petsc
  # always do this last
  setPCCtx(pc, mesh, sbp, eqn, opts, ctx_residual, t)

  return nothing
end

function applyPC(pc::PetscMatFreePC, mesh::AbstractMesh, sbp::AbstractSBP,
                 eqn::AbstractSolutionData, opts::Dict, t, b::AbstractVector, 
                 x::AbstractVector)

  # use the data in pc to multiply the preconditioner by b, storing result in x

  return nothing
end

function applyPCTranspose(pc::PetscMatFreePC, mesh::AbstractMesh,
                 sbp::AbstractSBP,
                 eqn::AbstractSolutionData, opts::Dict, t, b::AbstractVector, 
                 x::AbstractVector)

  # use the data in pc to multiply the transpose of the preconditioner by b,
  # storing result in x

  return nothing
end

Matrix-free preconditioners support the same kind of composition as matrix explicit preconditioners, so it is recommended to add the result into x, rather than overwrite x when pc_inner is not one of the Concrete PC Types.