17 #ifndef __deal2__matrix_tools_h 18 #define __deal2__matrix_tools_h 21 #include <deal.II/base/config.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> 36 template<
typename number>
class Vector;
43 template <
int dim,
int spacedim>
class Mapping;
44 template <
int dim,
int spacedim>
class DoFHandler;
46 template <
int dim,
int spacedim>
class FEValues;
56 #ifdef DEAL_II_WITH_PETSC 71 #ifdef DEAL_II_WITH_TRILINOS 242 template <
int dim,
typename number,
int spacedim>
254 template <
int dim,
typename number,
int spacedim>
274 template <
int dim,
typename number,
int spacedim>
288 template <
int dim,
typename number,
int spacedim>
300 template <
int dim,
typename number,
int spacedim>
311 template <
int dim,
typename number,
int spacedim>
321 template <
int dim,
typename number,
int spacedim>
334 template <
int dim,
typename number,
int spacedim>
366 template <
int dim,
int spacedim>
374 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
376 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
406 template <
int dim,
int spacedim>
413 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
415 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
420 template <
int dim,
int spacedim>
428 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
430 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
448 template <
int dim,
int spacedim>
455 std::vector<types::global_dof_index> &dof_to_boundary_mapping,
457 std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
472 template <
int dim,
int spacedim>
484 template <
int dim,
int spacedim>
504 template <
int dim,
int spacedim>
518 template <
int dim,
int spacedim>
531 template <
int dim,
int spacedim>
543 template <
int dim,
int spacedim>
554 template <
int dim,
int spacedim>
568 template <
int dim,
int spacedim>
804 template <
typename number>
810 const bool eliminate_columns =
true);
821 template <
typename number>
823 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
827 const bool eliminate_columns =
true);
829 #ifdef DEAL_II_WITH_PETSC 868 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
872 const bool eliminate_columns =
true);
879 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
883 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);
921 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
925 const bool eliminate_columns =
true);
929 #ifdef DEAL_II_WITH_TRILINOS 969 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
973 const bool eliminate_columns =
true);
981 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
985 const bool eliminate_columns =
true);
1024 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
1028 const bool eliminate_columns =
true);
1036 apply_boundary_values (
const std::map<types::global_dof_index,double> &boundary_values,
1040 const bool eliminate_columns =
true);
1081 const std::vector<types::global_dof_index> &local_dof_indices,
1084 const bool eliminate_columns);
1094 DEAL_II_NAMESPACE_CLOSE
std::map< types::boundary_id, const Function< dim > * > type
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())
#define DeclException0(Exception0)
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())