16 #ifndef dealii_matrix_tools_h 17 #define dealii_matrix_tools_h 20 #include <deal.II/base/config.h> 21 #include <deal.II/base/function.h> 22 #include <deal.II/base/exceptions.h> 23 #include <deal.II/base/thread_management.h> 24 #include <deal.II/lac/constraint_matrix.h> 25 #include <deal.II/dofs/function_map.h> 29 #ifdef DEAL_II_WITH_PETSC 30 # include <petscsys.h> 33 DEAL_II_NAMESPACE_OPEN
40 template <
typename number>
class Vector;
47 template <
int dim,
int spacedim>
class Mapping;
48 template <
int dim,
int spacedim>
class DoFHandler;
52 template <
int>
class QCollection;
53 template <
int,
int>
class MappingCollection;
58 #ifdef DEAL_II_WITH_PETSC 71 #ifdef DEAL_II_WITH_TRILINOS 238 template <
int dim,
int spacedim,
typename number>
250 template <
int dim,
int spacedim,
typename number>
276 template <
int dim,
int spacedim,
typename number>
290 template <
int dim,
int spacedim,
typename number>
302 template <
int dim,
int spacedim,
typename number>
313 template <
int dim,
int spacedim,
typename number>
323 template <
int dim,
int spacedim,
typename number>
336 template <
int dim,
int spacedim,
typename number>
371 template <
int dim,
int spacedim,
typename number>
378 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
380 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
387 template <
int dim,
int spacedim,
typename number>
393 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
395 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
400 template <
int dim,
int spacedim,
typename number>
407 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
409 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
414 template <
int dim,
int spacedim,
typename number>
420 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
422 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
443 template <
int dim,
int spacedim>
455 template <
int dim,
int spacedim>
480 template <
int dim,
int spacedim>
494 template <
int dim,
int spacedim>
507 template <
int dim,
int spacedim>
519 template <
int dim,
int spacedim>
530 template <
int dim,
int spacedim>
544 template <
int dim,
int spacedim>
557 "You are providing either a right hand side function or a " 558 "coefficient with a number of vector components that is " 559 "inconsistent with the rest of the arguments. If you do " 560 "provide a coefficient or right hand side function, then " 561 "it either needs to have as many components as the finite " 562 "element in use, or only a single vector component. In " 563 "the latter case, the same value will be taken for " 564 "each vector component of the finite element.");
773 template <
typename number>
779 const bool eliminate_columns =
true);
786 template <
typename number>
792 const bool eliminate_columns =
true);
794 #ifdef DEAL_II_WITH_PETSC 832 (
const std::map<types::global_dof_index,PetscScalar> &boundary_values,
836 const bool eliminate_columns =
true);
846 const bool eliminate_columns =
true);
850 #ifdef DEAL_II_WITH_TRILINOS 889 const bool eliminate_columns =
true);
900 const bool eliminate_columns =
true);
921 template <
typename number>
924 const std::vector<types::global_dof_index> &local_dof_indices,
927 const bool eliminate_columns);
933 "You are providing a matrix whose subdivision into " 934 "blocks in either row or column direction does not use " 935 "the same blocks sizes as the solution vector or " 936 "right hand side vectors, respectively.");
941 DEAL_II_NAMESPACE_CLOSE
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 ConstraintMatrix &constraints=ConstraintMatrix())
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 ConstraintMatrix &constraints=ConstraintMatrix())
Abstract base class for mapping classes.
#define DeclExceptionMsg(Exception, defaulttext)
static ::ExceptionBase & ExcComponentMismatch()
static ::ExceptionBase & ExcBlocksDontMatch()
void create_boundary_mass_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim-1 > &q, SparseMatrix< typename numbers::NumberTraits< number >::real_type > &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, typename numbers::NumberTraits< number >::real_type > *const weight=0, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())