16 #ifndef dealii_matrix_tools_h 17 #define dealii_matrix_tools_h 20 #include <deal.II/base/config.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/function.h> 24 #include <deal.II/base/thread_management.h> 26 #include <deal.II/dofs/deprecated_function_map.h> 28 #include <deal.II/lac/affine_constraints.h> 32 #ifdef DEAL_II_WITH_PETSC 33 # include <petscsys.h> 36 DEAL_II_NAMESPACE_OPEN
44 template <
typename number>
46 template <
typename number>
48 template <
typename number>
51 template <
typename number>
53 template <
typename Number>
56 template <
int dim,
int spacedim>
58 template <
int dim,
int spacedim>
66 class MappingCollection;
72 #ifdef DEAL_II_WITH_PETSC 85 #ifdef DEAL_II_WITH_TRILINOS 253 template <
int dim,
int spacedim,
typename number>
267 template <
int dim,
int spacedim,
typename number>
295 template <
int dim,
int spacedim,
typename number>
311 template <
int dim,
int spacedim,
typename number>
325 template <
int dim,
int spacedim,
typename number>
338 template <
int dim,
int spacedim,
typename number>
350 template <
int dim,
int spacedim,
typename number>
365 template <
int dim,
int spacedim,
typename number>
402 template <
int dim,
int spacedim,
typename number>
410 & boundary_functions,
412 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
414 std::vector<unsigned int> component_mapping = {});
421 template <
int dim,
int spacedim,
typename number>
428 & boundary_functions,
430 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
432 std::vector<unsigned int> component_mapping = {});
437 template <
int dim,
int spacedim,
typename number>
445 & boundary_functions,
447 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
449 std::vector<unsigned int> component_mapping = {});
454 template <
int dim,
int spacedim,
typename number>
461 & boundary_functions,
463 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
465 std::vector<unsigned int> component_mapping = {});
486 template <
int dim,
int spacedim>
500 template <
int dim,
int spacedim>
527 template <
int dim,
int spacedim>
543 template <
int dim,
int spacedim>
558 template <
int dim,
int spacedim>
572 template <
int dim,
int spacedim>
585 template <
int dim,
int spacedim>
601 template <
int dim,
int spacedim>
616 "You are providing either a right hand side function or a " 617 "coefficient with a number of vector components that is " 618 "inconsistent with the rest of the arguments. If you do " 619 "provide a coefficient or right hand side function, then " 620 "it either needs to have as many components as the finite " 621 "element in use, or only a single vector component. In " 622 "the latter case, the same value will be taken for " 623 "each vector component of the finite element.");
832 template <
typename number>
835 const std::map<types::global_dof_index, number> &boundary_values,
839 const bool eliminate_columns =
true);
846 template <
typename number>
849 const std::map<types::global_dof_index, number> &boundary_values,
853 const bool eliminate_columns =
true);
855 #ifdef DEAL_II_WITH_PETSC 893 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
897 const bool eliminate_columns =
true);
904 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
908 const bool eliminate_columns =
true);
912 #ifdef DEAL_II_WITH_TRILINOS 948 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
952 const bool eliminate_columns =
true);
960 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
964 const bool eliminate_columns =
true);
985 template <
typename number>
988 const std::map<types::global_dof_index, number> &boundary_values,
989 const std::vector<types::global_dof_index> & local_dof_indices,
992 const bool eliminate_columns);
998 "You are providing a matrix whose subdivision into " 999 "blocks in either row or column direction does not use " 1000 "the same blocks sizes as the solution vector or " 1001 "right hand side vectors, respectively.");
1006 DEAL_II_NAMESPACE_CLOSE
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_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={})
LinearAlgebra::distributed::Vector< Number > Vector
Abstract base class for mapping classes.
#define DeclExceptionMsg(Exception, defaulttext)
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 >())
static ::ExceptionBase & ExcComponentMismatch()
static ::ExceptionBase & ExcBlocksDontMatch()