DOLFINx 0.8.0
DOLFINx C++ interface
Loading...
Searching...
No Matches
Namespaces | Classes | Functions
dolfinx::la::petsc Namespace Reference

PETSc linear algebra functions. More...

Namespaces

namespace  options
 

Classes

class  KrylovSolver
 
class  Matrix
 
class  Operator
 
class  Vector
 

Functions

void error (int error_code, std::string filename, std::string petsc_function)
 Print error message for PETSc calls that return an error.
 
std::vector< Vec > create_vectors (MPI_Comm comm, const std::vector< std::span< const PetscScalar > > &x)
 
Vec create_vector (const common::IndexMap &map, int bs)
 
Vec create_vector (MPI_Comm comm, std::array< std::int64_t, 2 > range, std::span< const std::int64_t > ghosts, int bs)
 
Vec create_vector_wrap (const common::IndexMap &map, int bs, std::span< const PetscScalar > x)
 
template<class V >
Vec create_vector_wrap (const la::Vector< V > &x)
 
std::vector< IS > create_index_sets (const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
 
std::vector< std::vector< PetscScalar > > get_local_vectors (const Vec x, const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
 Copy blocks from Vec into local vectors.
 
void scatter_local_vectors (Vec x, const std::vector< std::span< const PetscScalar > > &x_b, const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > &maps)
 Scatter local vectors to Vec.
 
Mat create_matrix (MPI_Comm comm, const SparsityPattern &sp, std::string type=std::string())
 
MatNullSpace create_nullspace (MPI_Comm comm, std::span< const Vec > basis)
 

Detailed Description

PETSc linear algebra functions.

Function Documentation

◆ create_index_sets()

std::vector< IS > create_index_sets ( const std::vector< std::pair< std::reference_wrapper< const common::IndexMap >, int > > & maps)
Todo
This function could take just the local sizes

Compute PETSc IndexSets (IS) for a stack of index maps. E.g., if map[0] = {0, 1, 2, 3, 4, 5, 6} and map[1] = {0, 1, 2, 4} (in local indices) then IS[0] = {0, 1, 2, 3, 4, 5, 6} and IS[1] = {7, 8, 9, 10}.

Note
The caller is responsible for destruction of each IS.
Parameters
[in]mapsVector of IndexMaps and corresponding block sizes
Returns
Vector of PETSc Index Sets, created onPETSC_COMM_SELF

◆ create_matrix()

Mat create_matrix ( MPI_Comm comm,
const SparsityPattern & sp,
std::string type = std::string() )

Create a PETSc Mat. Caller is responsible for destroying the returned object.

◆ create_nullspace()

MatNullSpace create_nullspace ( MPI_Comm comm,
std::span< const Vec > basis )

Create PETSc MatNullSpace. Caller is responsible for destruction returned object.

Parameters
[in]commThe MPI communicator
[in]basisThe nullspace basis vectors
Returns
A PETSc nullspace object

◆ create_vector() [1/2]

Vec create_vector ( const common::IndexMap & map,
int bs )

Create a ghosted PETSc Vec

Note
Caller is responsible for destroying the returned object
Parameters
[in]mapThe index map describing the parallel layout (by block)
[in]bsThe block size
Returns
A PETSc Vec

◆ create_vector() [2/2]

Vec create_vector ( MPI_Comm comm,
std::array< std::int64_t, 2 > range,
std::span< const std::int64_t > ghosts,
int bs )

Create a ghosted PETSc Vec from a local range and ghost indices

Note
Caller is responsible for freeing the returned object
Parameters
[in]commThe MPI communicator
[in]rangeThe local ownership range (by blocks)
[in]ghostsGhost blocks
[in]bsThe block size. The total number of local entries is bs * (range[1] - range[0]).
Returns
A PETSc Vec

◆ create_vector_wrap() [1/2]

Vec create_vector_wrap ( const common::IndexMap & map,
int bs,
std::span< const PetscScalar > x )

Create a PETSc Vec that wraps the data in an array

Parameters
[in]mapThe index map that describes the parallel layout of the distributed vector (by block)
[in]bsBlock size
[in]xThe local part of the vector, including ghost entries
Returns
A PETSc Vec object that shares the data in x
Note
The array x must be kept alive to use the PETSc Vec object
The caller should call VecDestroy to free the return PETSc vector

◆ create_vector_wrap() [2/2]

template<class V >
Vec create_vector_wrap ( const la::Vector< V > & x)

Create a PETSc Vec that wraps the data in an array

Parameters
[in]xThe vector to be wrapped
Returns
A PETSc Vec object that shares the data in x

◆ create_vectors()

std::vector< Vec > create_vectors ( MPI_Comm comm,
const std::vector< std::span< const PetscScalar > > & x )

Create PETsc vectors from the local data. The data is copied into the PETSc vectors and is not shared.

Note
Caller is responsible for destroying the returned object
Parameters
[in]commThe MPI communicator
[in]xThe vector data owned by the calling rank. All components must have the same length.
Returns
Array of PETSc vectors