DOLFINx 0.9.0
DOLFINx C++ interface
|
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()) |
Create a matrix. | |
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. | |
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 alpha) |
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 alpha) |
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 alpha=1) |
Helper functions for assembly into PETSc data structures.
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 | alpha ) |
Modify b such that:
b <- b - alpha * 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.
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 | alpha ) |
Modify RHS vector to account for Dirichlet boundary conditions.
Modify b such that:
b <- b - alpha * 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.
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.
[in,out] | b | The 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] | L | The linear form to assemble |
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
.
[in,out] | b | The 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] | L | The linear form to assemble |
[in] | constants | The constants that appear in L |
[in] | coeffs | The coefficients that appear in L |
Mat create_matrix | ( | const Form< PetscScalar, T > & | a, |
std::string | type = std::string() ) |
Create a matrix.
[in] | a | A bilinear form |
[in] | type | The PETSc matrix type to create |
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.
[in] | a | Rectangular array of bilinear forms. The a(i, j) form will correspond to the (i, j) block in the returned matrix |
[in] | type | The type of PETSc Mat. If empty the PETSc default is used. |
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 | ) |
Initialise monolithic vector. Vector is not zeroed.
The caller is responsible for destroying the Mat object
void set_bc | ( | Vec | b, |
const std::vector< std::shared_ptr< const DirichletBC< PetscScalar, T > > > & | bcs, | ||
const Vec | x0, | ||
PetscScalar | alpha = 1 ) |
Set bc values in owned (local) part of the PETSc vector, multiplied by 'alpha'. 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.