fenapack package

Submodules

fenapack.assembling module

This module provides wrappers for assembling systems of linear algebraic equations to be solved with the use of PCD preconditioning strategy. These wrappers naturally provide routines for assembling preconditioning operators themselves.

class fenapack.assembling.PCDAssembler(a, L, bcs, a_pc=None, mp=None, mu=None, ap=None, fp=None, kp=None, gp=None, bcs_pcd=[])[source]

Bases: object

Base class for creating linear problems to be solved by application of the PCD preconditioning strategy. Users are encouraged to use this class for interfacing with fenapack.field_split.PCDKrylovSolver. On request it assembles not only the individual PCD operators but also the system matrix and the right hand side vector defining the linear problem.

__init__(a, L, bcs, a_pc=None, mp=None, mu=None, ap=None, fp=None, kp=None, gp=None, bcs_pcd=[])[source]

Collect individual variational forms and boundary conditions defining a linear problem (system matrix + RHS vector) on the one side and preconditioning operators on the other side.

Arguments
a (dolfin.Form or ufl.Form)
Bilinear form representing a system matrix.
L (dolfin.Form or ufl.Form)
Linear form representing a right hand side vector.
bcs (list of dolfin.DirichletBC)
Boundary conditions applied to a, L, and a_pc.
a_pc (dolfin.Form or ufl.Form)
Bilinear form representing a matrix optionally passed to preconditioner instead of a. In case of PCD, stabilized 00-block can be passed to 00-KSP solver.
mp, mu, ap, fp, kp, gp (dolfin.Form or ufl.Form)
Bilinear forms which (some of them) might be used by a particular PCD(R) preconditioner. Typically they represent “mass matrix” on pressure, “mass matrix” on velocity, minus Laplacian operator on pressure, pressure convection-diffusion operator, pressure convection operator and pressure gradient respectively.
bcs_pcd (list of dolfin.DirichletBC)
Artificial boundary conditions used by PCD preconditioner.

All the arguments should be given on the common mixed function space.

All the forms are wrapped using PCDForm so that each of them can be endowed with additional set of properties.

By default, mp, mu, ap and gp are assumed to be constant if the preconditioner is used repeatedly in some outer iterative process (e.g Newton-Raphson method, time-stepping). As such, the corresponding operators are assembled only once. On the other hand, fp and kp are updated in every outer iteration.

Also note that gp is the only form that is by default in a phantom mode. It means that the corresponding operator (if needed) is not obtained by assembling the form, but it is extracted as the 01-block of the system matrix.

The default setting can be modified by accessing a PCDForm instance via PCDAssembler.get_pcd_form() and changing the properties directly.

__weakref__

list of weak references to the object (if defined)

get_dolfin_form(key)[source]

Return form as dolfin.Form or ufl.Form.

get_pcd_form(key)[source]

Return form wrapped in PCDForm.

gp(Bt)[source]

Assemble discrete pressure gradient. It is crucial to respect any constraints placed on the velocity test space by Dirichlet boundary conditions.

pc_matrix(P)[source]

Assemble preconditioning matrix P whose relevant blocks can be passed to actual parts of the KSP solver.

rhs_vector(b, x=None)[source]

Assemble right hand side vector b.

The version with x is suitable for use inside a (quasi)-Newton solver.

system_matrix(A)[source]

Assemble system matrix A.

class fenapack.assembling.PCDForm(form, const=False, phantom=False)[source]

Bases: object

Wrapper for PCD operators represented by dolfin.Form or ufl.Form. This class allows to record specific properties of the form that can be utilized later while setting up the preconditioner.

For example, we can specify which matrices remain constant during the outer iterative algorithm (e.g. Newton-Raphson method, time-stepping) and which matrices need to be updated in every outer iteration.

__init__(form, const=False, phantom=False)[source]

The class is initialized by a single form with default properties.

Arguments
form (dolfin.Form or ufl.Form)
A form to be wrapped.
const (bool)
Whether the form remains constant in outer iterations.
phantom (bool)
If True, then the corresponding operator will be obtained not by assembling the form, but in some different way. (For example, pressure gradient may be extracted directly from the system matrix.)
__weakref__

list of weak references to the object (if defined)

fenapack.field_split module

This module provides subclasses of DOLFIN and petsc4py Krylov solvers implementing PCD fieldsplit preconditioned GMRES

fenapack.field_split_backend module

Tools for extraction and management of fieldsplit submatrices, subvectors, subbcs, subksps intended to be hidden from user interface

class fenapack.field_split_backend.PCDInterface(pcd_assembler, A, is_u, is_p, deep_submats=False)[source]

Bases: object

Wrapper of PCDAssembler for interfacing with PCD PC fieldsplit implementation. Convection fieldsplit submatrices are extracted as shallow or deep submatrices according to deep_submats parameter.

__init__(pcd_assembler, A, is_u, is_p, deep_submats=False)[source]

Create PCDInterface instance given PCDAssembler instance, system matrix and velocity and pressure index sets

__weakref__

list of weak references to the object (if defined)

apply_bcs(vec, bcs_getter, iset)[source]

Transform dolfin bcs obtained using bcs_getter function into fieldsplit subBCs and apply them to fieldsplit vector. SubBCs are cached.

apply_pcd_bcs(vec)[source]

Apply bcs to intermediate pressure vector of PCD pc

get_work_dolfin_mat(key, comm, can_be_destroyed=None, can_be_shared=None)[source]

Get working DOLFIN matrix by key. can_be_destroyed=True tells that it is probably favourable to not store the matrix unless it is shared as it will not be used ever again, None means that it can be destroyed but it is not probably favourable and False forbids the destruction. can_be_shared tells if a work matrix can be the same with work matrices for other keys.

get_work_mats(M, num)[source]

Return num of work mats initially created from matrix B.

get_work_vecs_from_square_mat(M, num)[source]

Return num of work vecs initially created from a square matrix M.

setup_ksp(ksp, assemble_func, iset, spd=False, const=False)[source]

Assemble into operator of given ksp if not yet assembled

setup_ksp_Ap(ksp)[source]

Setup pressure Laplacian ksp and assemble matrix

setup_ksp_Mp(ksp)[source]

Setup pressure mass matrix ksp and assemble matrix

setup_ksp_Rp(ksp, Mu, Bt)[source]

Setup pressure Laplacian ksp based on velocity mass matrix Mu and discrete gradient Bt and assemble matrix

setup_mat_Bt(mat=None)[source]

Setup and assemble discrete pressure gradient and return it

setup_mat_Fp(mat=None)[source]

Setup and assemble pressure convection-diffusion matrix and return it

setup_mat_Kp(mat=None)[source]

Setup and assemble pressure convection matrix and return it

setup_mat_Mu(mat=None)[source]

Setup and assemble velocity mass matrix and return it

fenapack.nonlinear_solvers module

This module provides subclasses of DOLFIN interface for solving non-linear problems suitable for use with fieldsplit preconditioned Krylov methods

fenapack.preconditioners module

class fenapack.preconditioners.BasePCDPC[source]

Bases: object

Base python context for pressure convection diffusion (PCD) preconditioners.

__weakref__

list of weak references to the object (if defined)

static create_default_ksp(comm)[source]

Return Cholesky factorization KSP

get_work_vecs(v, num)[source]

Return num of work vecs initially duplicated from v

init_pcd(pcd_interface)[source]

Initialize by PCDInterface instance

class fenapack.preconditioners.BasePCDRPC[source]

Bases: fenapack.preconditioners.BasePCDPC

Base python context for pressure convection diffusion reaction (PCDR) preconditioners.

class fenapack.preconditioners.PCDPC_BRM1[source]

Bases: fenapack.preconditioners.BasePCDPC

This class implements a modification of PCD variant similar to one by [1].

[1]Olshanskii M. A., Vassilevski Y. V., Pressure Schur complement preconditioners for the discrete Oseen problem. SIAM J. Sci. Comput., 29(6), 2686-2704. 2007.
class fenapack.preconditioners.PCDPC_BRM2[source]

Bases: fenapack.preconditioners.BasePCDPC

This class implements a modification of steady variant of PCD described in [2].

[2]Elman H. C., Silvester D. J., Wathen A. J., Finite Elements and Fast Iterative Solvers: With Application in Incompressible Fluid Dynamics. Oxford University Press 2005. 2nd edition 2014.
class fenapack.preconditioners.PCDRPC_BRM1[source]

Bases: fenapack.preconditioners.BasePCDRPC

This class implements an extension of PCDPC_BRM1. Here we add a reaction term into the preconditioner, so that is becomes PCDR (pressure-convection-diffusion-reaction) preconditioner. This particular variant is suitable for time-dependent problems, where the reaction term arises from the time derivative in the balance of momentum.

class fenapack.preconditioners.PCDRPC_BRM2[source]

Bases: fenapack.preconditioners.BasePCDRPC

This class implements an extension of PCDPC_BRM2. Here we add a reaction term into the preconditioner, so that is becomes PCDR (pressure-convection-diffusion-reaction) preconditioner. This particular variant is suitable for time-dependent problems, where the reaction term arises from the time derivative in the balance of momentum.

fenapack.stabilization module

fenapack.stabilization.StabilizationParameterSD(wind, viscosity, density=None)[source]

Returns a subclass of dolfin.Expression representing streamline diffusion stabilization parameter.

This kind of stabilization is convenient when a multigrid method is used for the convection term in the Navier-Stokes equation. The idea of the stabilization involves adding an additional term of the form:

delta_sd*inner(dot(grad(u), w), dot(grad(v), w))*dx

into the Navier-Stokes equation. Here u is a trial function, v is a test function and w defines so-called “wind” which is a known vector function. Regularization parameter delta_sd is determined by the local mesh Peclet number (PE), see the implementation below.

Arguments
wind (dolfin.GenericFunction)
A vector field determining convective velocity.
viscosity (dolfin.GenericFunction)
A scalar field determining dynamic viscosity.
density (dolfin.GenericFunction)
A scalar field determining density (optional).

fenapack.utils module

fenapack.utils.allow_only_one_call(func)[source]

Decorator allowing provided instancemethod to be called only once. Additional calls raise error.

fenapack.utils.get_default_factor_solver_type(comm)[source]

Return first available factor solver type name. This is implemened using DOLFIN now.

Module contents

This is FENaPack, FEniCS Navier-Stokes preconditioning package.