DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Functions
dolfinx::fem::petsc Namespace Reference

Helper functions for assembly into PETSc data structures. More...

Functions

template<std::floating_point T>
Mat create_matrix (const Form< PetscScalar, T > &a, std::string type=std::string())
 
template<std::floating_point T>
Mat create_matrix_block (const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, std::string type=std::string())
 
template<std::floating_point T>
Mat create_matrix_nest (const std::vector< std::vector< const Form< PetscScalar, T > * > > &a, const std::vector< std::vector< std::string > > &types)
 Create nested (MatNest) matrix.
 
Vec create_vector_block (const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
 
Vec create_vector_nest (const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
 Create nested (VecNest) vector. Vector is not zeroed.
 
template<std::floating_point T>
void assemble_vector (Vec b, const Form< PetscScalar, T > &L, std::span< const PetscScalar > constants, const std::map< std::pair< IntegralType, int >, std::pair< std::span< const PetscScalar >, int > > &coeffs)
 Assemble linear form into an already allocated PETSc vector.
 
template<std::floating_point T>
void assemble_vector (Vec b, const Form< PetscScalar, T > &L)
 
template<std::floating_point T>
void apply_lifting (Vec b, const std::vector< std::shared_ptr< const Form< PetscScalar, T > > > &a, const std::vector< std::span< const PetscScalar > > &constants, const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const PetscScalar >, int > > > &coeffs, const std::vector< std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > > &bcs1, const std::vector< Vec > &x0, PetscScalar scale)
 Modify RHS vector to account for Dirichlet boundary conditions.
 
template<std::floating_point T>
void apply_lifting (Vec b, const std::vector< std::shared_ptr< const Form< PetscScalar, double > > > &a, const std::vector< std::vector< std::shared_ptr< const DirichletBC< PetscScalar, double > > > > &bcs1, const std::vector< Vec > &x0, PetscScalar scale)
 
template<std::floating_point T>
void set_bc (Vec b, const std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > &bcs, const Vec x0, PetscScalar scale=1)
 

Detailed Description

Helper functions for assembly into PETSc data structures.

Function Documentation

◆ apply_lifting() [1/2]

template<std::floating_point T>
void apply_lifting ( Vec b,
const std::vector< std::shared_ptr< const Form< PetscScalar, double > > > & a,
const std::vector< std::vector< std::shared_ptr< const DirichletBC< PetscScalar, double > > > > & bcs1,
const std::vector< Vec > & x0,
PetscScalar scale )

Modify b such that:

b <- b - scale * A_j (g_j - x0_j)

where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs1 are on the trial spaces V_j. The forms in [a] must have the same test space as L (from which b was built), but the trial space may differ. If x0 is not supplied, then it is treated as zero.

Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.

◆ apply_lifting() [2/2]

template<std::floating_point T>
void apply_lifting ( Vec b,
const std::vector< std::shared_ptr< const Form< PetscScalar, T > > > & a,
const std::vector< std::span< const PetscScalar > > & constants,
const std::vector< std::map< std::pair< IntegralType, int >, std::pair< std::span< const PetscScalar >, int > > > & coeffs,
const std::vector< std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > > & bcs1,
const std::vector< Vec > & x0,
PetscScalar scale )

Modify RHS vector to account for Dirichlet boundary conditions.

Modify b such that:

b <- b - scale * A_j (g_j - x0_j)

where j is a block (nest) index. For a non-blocked problem j = 0. The boundary conditions bcs1 are on the trial spaces V_j. The forms in [a] must have the same test space as L (from which b was built), but the trial space may differ. If x0 is not supplied, then it is treated as zero.

Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.

◆ assemble_vector() [1/2]

template<std::floating_point T>
void assemble_vector ( Vec b,
const Form< PetscScalar, T > & L )

Assemble linear form into an already allocated PETSc vector. Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.

Parameters
[in,out]bThe PETsc vector to assemble the form into. The vector must already be initialised with the correct size. The process-local contribution of the form is assembled into this vector. It is not zeroed before assembly.
[in]LThe linear form to assemble

◆ assemble_vector() [2/2]

template<std::floating_point T>
void assemble_vector ( Vec b,
const Form< PetscScalar, T > & L,
std::span< const PetscScalar > constants,
const std::map< std::pair< IntegralType, int >, std::pair< std::span< const PetscScalar >, int > > & coeffs )

Assemble linear form into an already allocated PETSc vector.

Ghost contributions are not accumulated (not sent to owner). Caller is responsible for calling VecGhostUpdateBegin/End.

Parameters
[in,out]bThe PETsc vector to assemble the form into. The vector must already be initialised with the correct size. The process-local contribution of the form is assembled into this vector. It is not zeroed before assembly.
[in]LThe linear form to assemble
[in]constantsThe constants that appear in L
[in]coeffsThe coefficients that appear in L

◆ create_matrix()

template<std::floating_point T>
Mat create_matrix ( const Form< PetscScalar, T > & a,
std::string type = std::string() )

Create a matrix

Parameters
[in]aA bilinear form
[in]typeThe PETSc matrix type to create
Returns
A sparse matrix with a layout and sparsity that matches the bilinear form. The caller is responsible for destroying the Mat object.

◆ create_matrix_block()

template<std::floating_point T>
Mat create_matrix_block ( const std::vector< std::vector< const Form< PetscScalar, T > * > > & a,
std::string type = std::string() )

Initialise a monolithic matrix for an array of bilinear forms

Parameters
[in]aRectangular array of bilinear forms. The a(i, j) form will correspond to the (i, j) block in the returned matrix
[in]typeThe type of PETSc Mat. If empty the PETSc default is used.
Returns
A sparse matrix with a layout and sparsity that matches the bilinear forms. The caller is responsible for destroying the Mat object.

◆ create_matrix_nest()

template<std::floating_point T>
Mat create_matrix_nest ( const std::vector< std::vector< const Form< PetscScalar, T > * > > & a,
const std::vector< std::vector< std::string > > & types )

Create nested (MatNest) matrix.

Note
The caller is responsible for destroying the Mat object.

◆ create_vector_block()

Vec create_vector_block ( const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > & maps)

Initialise monolithic vector. Vector is not zeroed.

The caller is responsible for destroying the Mat object

◆ set_bc()

template<std::floating_point T>
void set_bc ( Vec b,
const std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > & bcs,
const Vec x0,
PetscScalar scale = 1 )

Set bc values in owned (local) part of the PETSc vector, multiplied by 'scale'. The vectors b and x0 must have the same local size. The bcs should be on (sub-)spaces of the form L that b represents.