Reference documentation for deal.II version 8.4.2
matrix_tools.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2016 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii__matrix_tools_h
17 #define dealii__matrix_tools_h
18 
19 
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>
26 
27 #include <map>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 
32 // forward declarations
33 template <int dim> class Quadrature;
34 
35 
36 template<typename number> class Vector;
37 template<typename number> class FullMatrix;
38 template<typename number> class SparseMatrix;
39 
40 template <typename number> class BlockSparseMatrix;
41 template <typename Number> class BlockVector;
42 
43 template <int dim, int spacedim> class Mapping;
44 template <int dim, int spacedim> class DoFHandler;
45 template <int dim, int spacedim> class FEValues;
46 
47 namespace hp
48 {
49  template <int> class QCollection;
50  template <int, int> class MappingCollection;
51  template <int, int> class DoFHandler;
52 }
53 
54 
55 #ifdef DEAL_II_WITH_PETSC
56 namespace PETScWrappers
57 {
58  class SparseMatrix;
59  class Vector;
60  namespace MPI
61  {
62  class SparseMatrix;
63  class BlockSparseMatrix;
64  class Vector;
65  class BlockVector;
66  }
67 }
68 #endif
69 
70 #ifdef DEAL_II_WITH_TRILINOS
71 namespace TrilinosWrappers
72 {
73  class SparseMatrix;
74  class Vector;
75  class BlockSparseMatrix;
76  class BlockVector;
77  namespace MPI
78  {
79  class Vector;
80  class BlockVector;
81  }
82 }
83 #endif
84 
85 
217 namespace MatrixCreator
218 {
235  template <int dim, typename number, int spacedim>
236  void create_mass_matrix (const Mapping<dim, spacedim> &mapping,
237  const DoFHandler<dim,spacedim> &dof,
238  const Quadrature<dim> &q,
239  SparseMatrix<number> &matrix,
240  const Function<spacedim> *const a = 0,
241  const ConstraintMatrix &constraints = ConstraintMatrix());
242 
247  template <int dim, typename number, int spacedim>
248  void create_mass_matrix (const DoFHandler<dim,spacedim> &dof,
249  const Quadrature<dim> &q,
250  SparseMatrix<number> &matrix,
251  const Function<spacedim> *const a = 0,
252  const ConstraintMatrix &constraints = ConstraintMatrix());
253 
270  template <int dim, typename number, int spacedim>
271  void create_mass_matrix (const Mapping<dim, spacedim> &mapping,
272  const DoFHandler<dim,spacedim> &dof,
273  const Quadrature<dim> &q,
274  SparseMatrix<number> &matrix,
275  const Function<spacedim> &rhs,
276  Vector<number> &rhs_vector,
277  const Function<spacedim> *const a = 0,
278  const ConstraintMatrix &constraints = ConstraintMatrix());
279 
284  template <int dim, typename number, int spacedim>
285  void create_mass_matrix (const DoFHandler<dim,spacedim> &dof,
286  const Quadrature<dim> &q,
287  SparseMatrix<number> &matrix,
288  const Function<spacedim> &rhs,
289  Vector<number> &rhs_vector,
290  const Function<spacedim> *const a = 0,
291  const ConstraintMatrix &constraints = ConstraintMatrix());
292 
296  template <int dim, typename number, int spacedim>
297  void create_mass_matrix (const hp::MappingCollection<dim,spacedim> &mapping,
298  const hp::DoFHandler<dim,spacedim> &dof,
299  const hp::QCollection<dim> &q,
300  SparseMatrix<number> &matrix,
301  const Function<spacedim> *const a = 0,
302  const ConstraintMatrix &constraints = ConstraintMatrix());
303 
307  template <int dim, typename number, int spacedim>
308  void create_mass_matrix (const hp::DoFHandler<dim,spacedim> &dof,
309  const hp::QCollection<dim> &q,
310  SparseMatrix<number> &matrix,
311  const Function<spacedim> *const a = 0,
312  const ConstraintMatrix &constraints = ConstraintMatrix());
313 
317  template <int dim, typename number, int spacedim>
318  void create_mass_matrix (const hp::MappingCollection<dim,spacedim> &mapping,
319  const hp::DoFHandler<dim,spacedim> &dof,
320  const hp::QCollection<dim> &q,
321  SparseMatrix<number> &matrix,
322  const Function<spacedim> &rhs,
323  Vector<number> &rhs_vector,
324  const Function<spacedim> *const a = 0,
325  const ConstraintMatrix &constraints = ConstraintMatrix());
326 
330  template <int dim, typename number, int spacedim>
331  void create_mass_matrix (const hp::DoFHandler<dim,spacedim> &dof,
332  const hp::QCollection<dim> &q,
333  SparseMatrix<number> &matrix,
334  const Function<spacedim> &rhs,
335  Vector<number> &rhs_vector,
336  const Function<spacedim> *const a = 0,
337  const ConstraintMatrix &constraints = ConstraintMatrix());
338 
339 
362  template <int dim, int spacedim>
364  const DoFHandler<dim,spacedim> &dof,
365  const Quadrature<dim-1> &q,
366  SparseMatrix<double> &matrix,
367  const typename FunctionMap<spacedim>::type &boundary_functions,
368  Vector<double> &rhs_vector,
369  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
370  const Function<spacedim> *const weight = 0,
371  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
372 
373 
378  template <int dim, int spacedim>
379  void create_boundary_mass_matrix (const DoFHandler<dim,spacedim> &dof,
380  const Quadrature<dim-1> &q,
381  SparseMatrix<double> &matrix,
382  const typename FunctionMap<spacedim>::type &boundary_functions,
383  Vector<double> &rhs_vector,
384  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
385  const Function<spacedim> *const a = 0,
386  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
387 
391  template <int dim, int spacedim>
392  void create_boundary_mass_matrix (const hp::MappingCollection<dim,spacedim> &mapping,
393  const hp::DoFHandler<dim,spacedim> &dof,
394  const hp::QCollection<dim-1> &q,
395  SparseMatrix<double> &matrix,
396  const typename FunctionMap<spacedim>::type &boundary_functions,
397  Vector<double> &rhs_vector,
398  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
399  const Function<spacedim> *const a = 0,
400  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
401 
405  template <int dim, int spacedim>
406  void create_boundary_mass_matrix (const hp::DoFHandler<dim,spacedim> &dof,
407  const hp::QCollection<dim-1> &q,
408  SparseMatrix<double> &matrix,
409  const typename FunctionMap<spacedim>::type &boundary_functions,
410  Vector<double> &rhs_vector,
411  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
412  const Function<spacedim> *const a = 0,
413  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
414 
431  template <int dim, int spacedim>
432  void create_laplace_matrix (const Mapping<dim, spacedim> &mapping,
433  const DoFHandler<dim,spacedim> &dof,
434  const Quadrature<dim> &q,
435  SparseMatrix<double> &matrix,
436  const Function<spacedim> *const a = 0,
437  const ConstraintMatrix &constraints = ConstraintMatrix());
438 
443  template <int dim, int spacedim>
444  void create_laplace_matrix (const DoFHandler<dim,spacedim> &dof,
445  const Quadrature<dim> &q,
446  SparseMatrix<double> &matrix,
447  const Function<spacedim> *const a = 0,
448  const ConstraintMatrix &constraints = ConstraintMatrix());
449 
465  template <int dim, int spacedim>
466  void create_laplace_matrix (const Mapping<dim, spacedim> &mapping,
467  const DoFHandler<dim,spacedim> &dof,
468  const Quadrature<dim> &q,
469  SparseMatrix<double> &matrix,
470  const Function<spacedim> &rhs,
471  Vector<double> &rhs_vector,
472  const Function<spacedim> *const a = 0,
473  const ConstraintMatrix &constraints = ConstraintMatrix());
474 
479  template <int dim, int spacedim>
480  void create_laplace_matrix (const DoFHandler<dim,spacedim> &dof,
481  const Quadrature<dim> &q,
482  SparseMatrix<double> &matrix,
483  const Function<spacedim> &rhs,
484  Vector<double> &rhs_vector,
485  const Function<spacedim> *const a = 0,
486  const ConstraintMatrix &constraints = ConstraintMatrix());
487 
492  template <int dim, int spacedim>
493  void create_laplace_matrix (const hp::MappingCollection<dim,spacedim> &mapping,
494  const hp::DoFHandler<dim,spacedim> &dof,
495  const hp::QCollection<dim> &q,
496  SparseMatrix<double> &matrix,
497  const Function<spacedim> *const a = 0,
498  const ConstraintMatrix &constraints = ConstraintMatrix());
499 
504  template <int dim, int spacedim>
505  void create_laplace_matrix (const hp::DoFHandler<dim,spacedim> &dof,
506  const hp::QCollection<dim> &q,
507  SparseMatrix<double> &matrix,
508  const Function<spacedim> *const a = 0,
509  const ConstraintMatrix &constraints = ConstraintMatrix());
510 
515  template <int dim, int spacedim>
516  void create_laplace_matrix (const hp::MappingCollection<dim,spacedim> &mapping,
517  const hp::DoFHandler<dim,spacedim> &dof,
518  const hp::QCollection<dim> &q,
519  SparseMatrix<double> &matrix,
520  const Function<spacedim> &rhs,
521  Vector<double> &rhs_vector,
522  const Function<spacedim> *const a = 0,
523  const ConstraintMatrix &constraints = ConstraintMatrix());
524 
529  template <int dim, int spacedim>
530  void create_laplace_matrix (const hp::DoFHandler<dim,spacedim> &dof,
531  const hp::QCollection<dim> &q,
532  SparseMatrix<double> &matrix,
533  const Function<spacedim> &rhs,
534  Vector<double> &rhs_vector,
535  const Function<spacedim> *const a = 0,
536  const ConstraintMatrix &constraints = ConstraintMatrix());
537 
541  DeclExceptionMsg (ExcComponentMismatch,
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.");
550 }
551 
552 
553 
745 namespace MatrixTools
746 {
752  using namespace MatrixCreator;
753 
758  template <typename number>
759  void
760  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
761  SparseMatrix<number> &matrix,
762  Vector<number> &solution,
763  Vector<number> &right_hand_side,
764  const bool eliminate_columns = true);
765 
771  template <typename number>
772  void
773  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
775  BlockVector<number> &solution,
776  BlockVector<number> &right_hand_side,
777  const bool eliminate_columns = true);
778 
779 #ifdef DEAL_II_WITH_PETSC
780 
799  void
800  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
802  PETScWrappers::Vector &solution,
803  PETScWrappers::Vector &right_hand_side,
804  const bool eliminate_columns = true);
805 
825  void
826  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
828  PETScWrappers::MPI::Vector &solution,
829  PETScWrappers::MPI::Vector &right_hand_side,
830  const bool eliminate_columns = true);
831 
835  void
836  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
839  PETScWrappers::MPI::BlockVector &right_hand_side,
840  const bool eliminate_columns = true);
841 
842 #endif
843 
844 #ifdef DEAL_II_WITH_TRILINOS
845 
862  void
863  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
865  TrilinosWrappers::Vector &solution,
866  TrilinosWrappers::Vector &right_hand_side,
867  const bool eliminate_columns = true);
868 
873  void
874  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
877  TrilinosWrappers::BlockVector &right_hand_side,
878  const bool eliminate_columns = true);
879 
899  void
900  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
903  TrilinosWrappers::MPI::Vector &right_hand_side,
904  const bool eliminate_columns = true);
905 
910  void
911  apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
914  TrilinosWrappers::MPI::BlockVector &right_hand_side,
915  const bool eliminate_columns = true);
916 #endif
917 
936  void
937  local_apply_boundary_values (const std::map<types::global_dof_index,double> &boundary_values,
938  const std::vector<types::global_dof_index> &local_dof_indices,
939  FullMatrix<double> &local_matrix,
940  Vector<double> &local_rhs,
941  const bool eliminate_columns);
942 
946  DeclExceptionMsg (ExcBlocksDontMatch,
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.");
951 }
952 
953 
954 
955 DEAL_II_NAMESPACE_CLOSE
956 
957 #endif
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())
void local_apply_boundary_values(const std::map< types::global_dof_index, double > &boundary_values, const std::vector< types::global_dof_index > &local_dof_indices, FullMatrix< double > &local_matrix, Vector< double > &local_rhs, const bool eliminate_columns)
Abstract base class for mapping classes.
Definition: dof_tools.h:52
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:533
std::map< types::boundary_id, const Function< dim, Number > * > type
Definition: function_map.h:81
Definition: hp.h:102
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 >())
Definition: fe.h:31
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())
void apply_boundary_values(const std::map< types::global_dof_index, double > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
Definition: matrix_tools.cc:77