Inherits Base< dim, ::LinearAlgebra::distributed::BlockVector< number > >.
|
void | apply_add (::LinearAlgebra::distributed::BlockVector< number > &dst, const ::LinearAlgebra::distributed::BlockVector< number > &src) const override |
|
void | local_apply (const ::MatrixFree< dim, number > &data, ::LinearAlgebra::distributed::BlockVector< number > &dst, const ::LinearAlgebra::distributed::BlockVector< number > &src, const std::pair< unsigned int, unsigned int > &cell_range) const |
|
void | local_apply_face (const ::MatrixFree< dim, number > &data, ::LinearAlgebra::distributed::BlockVector< number > &dst, const ::LinearAlgebra::distributed::BlockVector< number > &src, const std::pair< unsigned int, unsigned int > &face_range) const |
|
void | local_apply_boundary_face (const ::MatrixFree< dim, number > &data, ::LinearAlgebra::distributed::BlockVector< number > &dst, const ::LinearAlgebra::distributed::BlockVector< number > &src, const std::pair< unsigned int, unsigned int > &face_range) const |
|
template<int dim, int degree_v, typename number>
class aspect::MatrixFreeStokesOperators::StokesOperator< dim, degree_v, number >
Operator for the entire Stokes block.
Definition at line 188 of file stokes_matrix_free.h.
§ StokesOperator()
template<int dim, int degree_v, typename number>
§ clear()
template<int dim, int degree_v, typename number>
§ set_cell_data()
template<int dim, int degree_v, typename number>
Pass in a reference to the problem data.
§ compute_diagonal()
template<int dim, int degree_v, typename number>
Computes the diagonal of the matrix. Since matrix-free operators have not access to matrix elements, we must apply the matrix-free operator to the unit vectors to recover the diagonal.
§ apply_add()
template<int dim, int degree_v, typename number>
Performs the application of the matrix-free operator. This function is called by vmult() functions MatrixFreeOperators::Base.
§ local_apply()
template<int dim, int degree_v, typename number>
void aspect::MatrixFreeStokesOperators::StokesOperator< dim, degree_v, number >::local_apply |
( |
const ::MatrixFree< dim, number > & |
data, |
|
|
::LinearAlgebra::distributed::BlockVector< number > & |
dst, |
|
|
const ::LinearAlgebra::distributed::BlockVector< number > & |
src, |
|
|
const std::pair< unsigned int, unsigned int > & |
cell_range |
|
) |
| const |
|
private |
Defines the application of the cell matrix.
§ local_apply_face()
template<int dim, int degree_v, typename number>
void aspect::MatrixFreeStokesOperators::StokesOperator< dim, degree_v, number >::local_apply_face |
( |
const ::MatrixFree< dim, number > & |
data, |
|
|
::LinearAlgebra::distributed::BlockVector< number > & |
dst, |
|
|
const ::LinearAlgebra::distributed::BlockVector< number > & |
src, |
|
|
const std::pair< unsigned int, unsigned int > & |
face_range |
|
) |
| const |
|
private |
This function doesn't do anything, it's created to use the matrixfree loop.
§ local_apply_boundary_face()
template<int dim, int degree_v, typename number>
void aspect::MatrixFreeStokesOperators::StokesOperator< dim, degree_v, number >::local_apply_boundary_face |
( |
const ::MatrixFree< dim, number > & |
data, |
|
|
::LinearAlgebra::distributed::BlockVector< number > & |
dst, |
|
|
const ::LinearAlgebra::distributed::BlockVector< number > & |
src, |
|
|
const std::pair< unsigned int, unsigned int > & |
face_range |
|
) |
| const |
|
private |
Apply the stabilization on free surface faces.
§ cell_data
template<int dim, int degree_v, typename number>
A pointer to the current cell data that contains viscosity and other required parameters per cell.
Definition at line 251 of file stokes_matrix_free.h.
The documentation for this class was generated from the following file: