16#ifndef dealii_matrix_tools_h
17#define dealii_matrix_tools_h
29#ifdef DEAL_II_WITH_PETSC
42template <
typename number>
44template <
typename number>
46template <
typename number>
49template <
typename number>
51template <
typename Number>
54template <
int dim,
int spacedim>
56template <
int dim,
int spacedim>
64 class MappingCollection;
70# ifdef DEAL_II_WITH_PETSC
83# ifdef DEAL_II_WITH_TRILINOS
251 template <
int dim,
int spacedim,
typename number>
265 template <
int dim,
int spacedim,
typename number>
293 template <
int dim,
int spacedim,
typename number>
309 template <
int dim,
int spacedim,
typename number>
323 template <
int dim,
int spacedim,
typename number>
336 template <
int dim,
int spacedim,
typename number>
348 template <
int dim,
int spacedim,
typename number>
363 template <
int dim,
int spacedim,
typename number>
400 template <
int dim,
int spacedim,
typename number>
408 & boundary_functions,
410 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
412 std::vector<unsigned int> component_mapping = {});
419 template <
int dim,
int spacedim,
typename number>
426 & boundary_functions,
428 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
430 std::vector<unsigned int> component_mapping = {});
435 template <
int dim,
int spacedim,
typename number>
443 & boundary_functions,
445 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
447 std::vector<unsigned int> component_mapping = {});
452 template <
int dim,
int spacedim,
typename number>
459 & boundary_functions,
461 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
463 std::vector<unsigned int> component_mapping = {});
484 template <
int dim,
int spacedim>
498 template <
int dim,
int spacedim>
525 template <
int dim,
int spacedim>
541 template <
int dim,
int spacedim>
555 template <
int dim,
int spacedim>
568 template <
int dim,
int spacedim>
580 template <
int dim,
int spacedim>
595 template <
int dim,
int spacedim>
610 "You are providing either a right hand side function or a "
611 "coefficient with a number of vector components that is "
612 "inconsistent with the rest of the arguments. If you do "
613 "provide a coefficient or right hand side function, then "
614 "it either needs to have as many components as the finite "
615 "element in use, or only a single vector component. In "
616 "the latter case, the same value will be taken for "
617 "each vector component of the finite element.");
842 template <
typename number>
845 const std::map<types::global_dof_index, number> &boundary_values,
849 const bool eliminate_columns =
true);
856 template <
typename number>
859 const std::map<types::global_dof_index, number> &boundary_values,
863 const bool eliminate_columns =
true);
865#ifdef DEAL_II_WITH_PETSC
903 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
907 const bool eliminate_columns =
true);
914 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
918 const bool eliminate_columns =
true);
922#ifdef DEAL_II_WITH_TRILINOS
958 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
962 const bool eliminate_columns =
true);
970 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
974 const bool eliminate_columns =
true);
995 template <
typename number>
998 const std::map<types::global_dof_index, number> &boundary_values,
999 const std::vector<types::global_dof_index> & local_dof_indices,
1002 const bool eliminate_columns);
1008 "You are providing a matrix whose subdivision into "
1009 "blocks in either row or column direction does not use "
1010 "the same blocks sizes as the solution vector or "
1011 "right hand side vectors, respectively.");
Abstract base class for mapping classes.
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcBlocksDontMatch()
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcComponentMismatch()
@ matrix
Contents is actually a matrix.
void create_boundary_mass_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim - 1 > &q, SparseMatrix< number > &matrix, const std::map< types::boundary_id, const Function< spacedim, number > * > &boundary_functions, Vector< number > &rhs_vector, std::vector< types::global_dof_index > &dof_to_boundary_mapping, const Function< spacedim, number > *const weight=0, std::vector< unsigned int > component_mapping={})
void create_laplace_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, SparseMatrix< double > &matrix, const Function< spacedim > *const a=nullptr, const AffineConstraints< double > &constraints=AffineConstraints< double >())
void create_mass_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, SparseMatrix< number > &matrix, const Function< spacedim, number > *const a=nullptr, const AffineConstraints< number > &constraints=AffineConstraints< number >())