Reference documentation for deal.II version 9.2.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
matrix_tools.h
Go to the documentation of this file.
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2020 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.md at
12 // the top level directory of deal.II.
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 
23 #include <deal.II/base/function.h>
24 
26 
28 
29 #include <map>
30 
31 #ifdef DEAL_II_WITH_PETSC
32 # include <petscsys.h>
33 #endif
34 
36 
37 
38 // forward declarations
39 #ifndef DOXYGEN
40 template <int dim>
41 class Quadrature;
42 
43 
44 template <typename number>
45 class Vector;
46 template <typename number>
47 class FullMatrix;
48 template <typename number>
49 class SparseMatrix;
50 
51 template <typename number>
52 class BlockSparseMatrix;
53 template <typename Number>
54 class BlockVector;
55 
56 template <int dim, int spacedim>
57 class Mapping;
58 template <int dim, int spacedim>
59 class DoFHandler;
60 
61 namespace hp
62 {
63  template <int>
64  class QCollection;
65  template <int, int>
66  class MappingCollection;
67  template <int, int>
68  class DoFHandler;
69 } // namespace hp
70 
71 
72 # ifdef DEAL_II_WITH_PETSC
73 namespace PETScWrappers
74 {
75  class MatrixBase;
76  class VectorBase;
77  namespace MPI
78  {
79  class BlockSparseMatrix;
80  class BlockVector;
81  } // namespace MPI
82 } // namespace PETScWrappers
83 # endif
84 
85 # ifdef DEAL_II_WITH_TRILINOS
86 namespace TrilinosWrappers
87 {
88  class SparseMatrix;
89  class BlockSparseMatrix;
90  namespace MPI
91  {
92  class Vector;
93  class BlockVector;
94  } // namespace MPI
95 } // namespace TrilinosWrappers
96 # endif
97 #endif
98 
99 
233 namespace MatrixCreator
234 {
254  template <int dim, int spacedim, typename number>
255  void
257  const Mapping<dim, spacedim> & mapping,
258  const DoFHandler<dim, spacedim> & dof,
259  const Quadrature<dim> & q,
261  const Function<spacedim, number> *const a = nullptr,
263 
268  template <int dim, int spacedim, typename number>
269  void
271  const DoFHandler<dim, spacedim> & dof,
272  const Quadrature<dim> & q,
274  const Function<spacedim, number> *const a = nullptr,
276 
296  template <int dim, int spacedim, typename number>
297  void
299  const Mapping<dim, spacedim> & mapping,
300  const DoFHandler<dim, spacedim> & dof,
301  const Quadrature<dim> & q,
303  const Function<spacedim, number> & rhs,
304  Vector<number> & rhs_vector,
305  const Function<spacedim, number> *const a = nullptr,
307 
312  template <int dim, int spacedim, typename number>
313  void
315  const DoFHandler<dim, spacedim> & dof,
316  const Quadrature<dim> & q,
318  const Function<spacedim, number> & rhs,
319  Vector<number> & rhs_vector,
320  const Function<spacedim, number> *const a = nullptr,
322 
326  template <int dim, int spacedim, typename number>
327  void
330  const hp::DoFHandler<dim, spacedim> & dof,
331  const hp::QCollection<dim> & q,
333  const Function<spacedim, number> *const a = nullptr,
335 
339  template <int dim, int spacedim, typename number>
340  void
342  const hp::DoFHandler<dim, spacedim> & dof,
343  const hp::QCollection<dim> & q,
345  const Function<spacedim, number> *const a = nullptr,
347 
351  template <int dim, int spacedim, typename number>
352  void
355  const hp::DoFHandler<dim, spacedim> & dof,
356  const hp::QCollection<dim> & q,
358  const Function<spacedim, number> & rhs,
359  Vector<number> & rhs_vector,
360  const Function<spacedim, number> *const a = nullptr,
362 
366  template <int dim, int spacedim, typename number>
367  void
369  const hp::DoFHandler<dim, spacedim> & dof,
370  const hp::QCollection<dim> & q,
372  const Function<spacedim, number> & rhs,
373  Vector<number> & rhs_vector,
374  const Function<spacedim, number> *const a = nullptr,
376 
377 
403  template <int dim, int spacedim, typename number>
404  void
406  const Mapping<dim, spacedim> & mapping,
407  const DoFHandler<dim, spacedim> &dof,
408  const Quadrature<dim - 1> & q,
410  const std::map<types::boundary_id, const Function<spacedim, number> *>
411  & boundary_functions,
412  Vector<number> & rhs_vector,
413  std::vector<types::global_dof_index> & dof_to_boundary_mapping,
414  const Function<spacedim, number> *const weight = 0,
415  std::vector<unsigned int> component_mapping = {});
416 
417 
422  template <int dim, int spacedim, typename number>
423  void
425  const DoFHandler<dim, spacedim> &dof,
426  const Quadrature<dim - 1> & q,
428  const std::map<types::boundary_id, const Function<spacedim, number> *>
429  & boundary_functions,
430  Vector<number> & rhs_vector,
431  std::vector<types::global_dof_index> & dof_to_boundary_mapping,
432  const Function<spacedim, number> *const a = nullptr,
433  std::vector<unsigned int> component_mapping = {});
434 
438  template <int dim, int spacedim, typename number>
439  void
442  const hp::DoFHandler<dim, spacedim> & dof,
443  const hp::QCollection<dim - 1> & q,
445  const std::map<types::boundary_id, const Function<spacedim, number> *>
446  & boundary_functions,
447  Vector<number> & rhs_vector,
448  std::vector<types::global_dof_index> & dof_to_boundary_mapping,
449  const Function<spacedim, number> *const a = nullptr,
450  std::vector<unsigned int> component_mapping = {});
451 
455  template <int dim, int spacedim, typename number>
456  void
459  const hp::QCollection<dim - 1> & q,
461  const std::map<types::boundary_id, const Function<spacedim, number> *>
462  & boundary_functions,
463  Vector<number> & rhs_vector,
464  std::vector<types::global_dof_index> & dof_to_boundary_mapping,
465  const Function<spacedim, number> *const a = nullptr,
466  std::vector<unsigned int> component_mapping = {});
467 
487  template <int dim, int spacedim>
488  void
490  const Mapping<dim, spacedim> & mapping,
491  const DoFHandler<dim, spacedim> &dof,
492  const Quadrature<dim> & q,
494  const Function<spacedim> *const a = nullptr,
496 
501  template <int dim, int spacedim>
502  void
504  const DoFHandler<dim, spacedim> &dof,
505  const Quadrature<dim> & q,
507  const Function<spacedim> *const a = nullptr,
509 
528  template <int dim, int spacedim>
529  void
531  const Mapping<dim, spacedim> & mapping,
532  const DoFHandler<dim, spacedim> &dof,
533  const Quadrature<dim> & q,
535  const Function<spacedim> & rhs,
536  Vector<double> & rhs_vector,
537  const Function<spacedim> *const a = nullptr,
539 
544  template <int dim, int spacedim>
545  void
547  const DoFHandler<dim, spacedim> &dof,
548  const Quadrature<dim> & q,
550  const Function<spacedim> & rhs,
551  Vector<double> & rhs_vector,
552  const Function<spacedim> *const a = nullptr,
554 
559  template <int dim, int spacedim>
560  void
563  const hp::DoFHandler<dim, spacedim> & dof,
564  const hp::QCollection<dim> & q,
566  const Function<spacedim> *const a = nullptr,
568 
573  template <int dim, int spacedim>
574  void
577  const hp::QCollection<dim> & q,
579  const Function<spacedim> *const a = nullptr,
581 
586  template <int dim, int spacedim>
587  void
590  const hp::DoFHandler<dim, spacedim> & dof,
591  const hp::QCollection<dim> & q,
593  const Function<spacedim> & rhs,
594  Vector<double> & rhs_vector,
595  const Function<spacedim> *const a = nullptr,
597 
602  template <int dim, int spacedim>
603  void
606  const hp::QCollection<dim> & q,
608  const Function<spacedim> & rhs,
609  Vector<double> & rhs_vector,
610  const Function<spacedim> *const a = nullptr,
612 
617  "You are providing either a right hand side function or a "
618  "coefficient with a number of vector components that is "
619  "inconsistent with the rest of the arguments. If you do "
620  "provide a coefficient or right hand side function, then "
621  "it either needs to have as many components as the finite "
622  "element in use, or only a single vector component. In "
623  "the latter case, the same value will be taken for "
624  "each vector component of the finite element.");
625 } // namespace MatrixCreator
626 
627 
628 
821 namespace MatrixTools
822 {
828  using namespace MatrixCreator;
829 
834  template <typename number>
835  void
837  const std::map<types::global_dof_index, number> &boundary_values,
839  Vector<number> & solution,
840  Vector<number> & right_hand_side,
841  const bool eliminate_columns = true);
842 
848  template <typename number>
849  void
851  const std::map<types::global_dof_index, number> &boundary_values,
853  BlockVector<number> & solution,
854  BlockVector<number> & right_hand_side,
855  const bool eliminate_columns = true);
856 
857 #ifdef DEAL_II_WITH_PETSC
858 
893  void
895  const std::map<types::global_dof_index, PetscScalar> &boundary_values,
897  PETScWrappers::VectorBase & solution,
898  PETScWrappers::VectorBase & right_hand_side,
899  const bool eliminate_columns = true);
900 
904  void
906  const std::map<types::global_dof_index, PetscScalar> &boundary_values,
909  PETScWrappers::MPI::BlockVector & right_hand_side,
910  const bool eliminate_columns = true);
911 
912 #endif
913 
914 #ifdef DEAL_II_WITH_TRILINOS
915 
948  void
950  const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
953  TrilinosWrappers::MPI::Vector & right_hand_side,
954  const bool eliminate_columns = true);
955 
960  void
962  const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
965  TrilinosWrappers::MPI::BlockVector & right_hand_side,
966  const bool eliminate_columns = true);
967 #endif
968 
987  template <typename number>
988  void
990  const std::map<types::global_dof_index, number> &boundary_values,
991  const std::vector<types::global_dof_index> & local_dof_indices,
992  FullMatrix<number> & local_matrix,
993  Vector<number> & local_rhs,
994  const bool eliminate_columns);
995 
1000  "You are providing a matrix whose subdivision into "
1001  "blocks in either row or column direction does not use "
1002  "the same blocks sizes as the solution vector or "
1003  "right hand side vectors, respectively.");
1004 } // namespace MatrixTools
1005 
1006 
1007 
1009 
1010 #endif
DeclExceptionMsg
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:496
TrilinosWrappers::MPI::Vector
Definition: trilinos_vector.h:400
LinearAlgebraDealII::Vector
Vector< double > Vector
Definition: generic_linear_algebra.h:43
deprecated_function_map.h
PETScWrappers
Definition: petsc_block_sparse_matrix.h:36
TrilinosWrappers::SparseMatrix
Definition: trilinos_sparse_matrix.h:515
BlockVector
Definition: block_vector.h:71
MatrixCreator::create_laplace_matrix
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 >())
SparseMatrix
Definition: sparse_matrix.h:497
TrilinosWrappers
Definition: types.h:161
hp::QCollection
Definition: q_collection.h:48
MatrixTools::apply_boundary_values
void apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, SparseMatrix< number > &matrix, Vector< number > &solution, Vector< number > &right_hand_side, const bool eliminate_columns=true)
Definition: matrix_tools.cc:81
TrilinosWrappers::MPI::BlockVector
Definition: trilinos_parallel_block_vector.h:75
hp::DoFHandler
Definition: dof_handler.h:203
DoFHandler
Definition: dof_handler.h:205
Mapping
Abstract base class for mapping classes.
Definition: mapping.h:302
MatrixCreator::create_boundary_mass_matrix
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={})
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
PETScWrappers::MatrixBase
Definition: petsc_matrix_base.h:284
MatrixCreator::create_mass_matrix
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 >())
LAPACKSupport::matrix
@ matrix
Contents is actually a matrix.
Definition: lapack_support.h:60
MatrixTools::ExcBlocksDontMatch
static ::ExceptionBase & ExcBlocksDontMatch()
hp::MappingCollection
Definition: mapping_collection.h:56
MatrixTools
Definition: matrix_tools.h:821
PETScWrappers::MPI::BlockVector
Definition: petsc_block_vector.h:61
exceptions.h
MatrixCreator
Definition: matrix_tools.h:233
unsigned int
AffineConstraints
Definition: affine_constraints.h:180
function.h
hp
Definition: hp.h:117
MatrixCreator::ExcComponentMismatch
static ::ExceptionBase & ExcComponentMismatch()
TrilinosWrappers::BlockSparseMatrix
Definition: trilinos_block_sparse_matrix.h:72
affine_constraints.h
config.h
FullMatrix
Definition: full_matrix.h:71
MatrixTools::local_apply_boundary_values
void local_apply_boundary_values(const std::map< types::global_dof_index, number > &boundary_values, const std::vector< types::global_dof_index > &local_dof_indices, FullMatrix< number > &local_matrix, Vector< number > &local_rhs, const bool eliminate_columns)
Definition: matrix_tools.cc:505
Function
Definition: function.h:151
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
Quadrature
Definition: quadrature.h:85
PETScWrappers::VectorBase
Definition: petsc_vector_base.h:240
Vector< number >
PETScWrappers::MPI::BlockSparseMatrix
Definition: petsc_block_sparse_matrix.h:67
BlockSparseMatrix
Definition: block_sparse_matrix.h:50