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 DEAL_II_NAMESPACE_OPEN
36 template<
typename number>
class Vector;
43 template <
int dim,
int spacedim>
class Mapping;
44 template <
int dim,
int spacedim>
class DoFHandler;
45 template <
int dim,
int spacedim>
class FEValues;
49 template <
int>
class QCollection;
50 template <
int,
int>
class MappingCollection;
55 #ifdef DEAL_II_WITH_PETSC 70 #ifdef DEAL_II_WITH_TRILINOS 235 template <
int dim,
typename number,
int spacedim>
247 template <
int dim,
typename number,
int spacedim>
270 template <
int dim,
typename number,
int spacedim>
284 template <
int dim,
typename number,
int spacedim>
296 template <
int dim,
typename number,
int spacedim>
307 template <
int dim,
typename number,
int spacedim>
317 template <
int dim,
typename number,
int spacedim>
330 template <
int dim,
typename number,
int spacedim>
362 template <
int dim,
int spacedim>
369 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
371 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
378 template <
int dim,
int spacedim>
384 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
386 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
391 template <
int dim,
int spacedim>
398 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
400 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
405 template <
int dim,
int spacedim>
411 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
413 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
431 template <
int dim,
int spacedim>
443 template <
int dim,
int spacedim>
465 template <
int dim,
int spacedim>
479 template <
int dim,
int spacedim>
492 template <
int dim,
int spacedim>
504 template <
int dim,
int spacedim>
515 template <
int dim,
int spacedim>
529 template <
int dim,
int spacedim>
542 "You are providing either a right hand side function or a " 543 "coefficient with a number of vector components that is " 544 "inconsistent with the rest of the arguments. If you do " 545 "provide a coefficient or right hand side function, then " 546 "it either needs to have as many components as the finite " 547 "element in use, or only a single vector component. In " 548 "the latter case, the same value will be taken for " 549 "each vector component of the finite element.");
758 template <
typename number>
764 const bool eliminate_columns =
true);
771 template <
typename number>
773 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
777 const bool eliminate_columns =
true);
779 #ifdef DEAL_II_WITH_PETSC 800 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
804 const bool eliminate_columns =
true);
826 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
830 const bool eliminate_columns =
true);
836 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
840 const bool eliminate_columns =
true);
844 #ifdef DEAL_II_WITH_TRILINOS 863 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
867 const bool eliminate_columns =
true);
874 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
878 const bool eliminate_columns =
true);
900 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
904 const bool eliminate_columns =
true);
911 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
915 const bool eliminate_columns =
true);
938 const std::vector<types::global_dof_index> &local_dof_indices,
941 const bool eliminate_columns);
947 "You are providing a matrix whose subdivision into " 948 "blocks in either row or column direction does not use " 949 "the same blocks sizes as the solution vector or " 950 "right hand side vectors, respectively.");
955 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=0, const ConstraintMatrix &constraints=ConstraintMatrix())
Abstract base class for mapping classes.
#define DeclExceptionMsg(Exception, defaulttext)
std::map< types::boundary_id, const Function< dim, Number > * > type
void create_boundary_mass_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim-1 > &q, SparseMatrix< double > &matrix, const typename FunctionMap< spacedim >::type &boundary_functions, Vector< double > &rhs_vector, std::vector< types::global_dof_index > &dof_to_boundary_mapping, const Function< spacedim > *const weight=0, std::vector< unsigned int > component_mapping=std::vector< unsigned int >())
void create_mass_matrix(const Mapping< dim, spacedim > &mapping, const DoFHandler< dim, spacedim > &dof, const Quadrature< dim > &q, SparseMatrix< number > &matrix, const Function< spacedim > *const a=0, const ConstraintMatrix &constraints=ConstraintMatrix())