Reference documentation for deal.II version 9.0.0
matrix_tools.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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 #ifdef DEAL_II_WITH_PETSC
30 # include <petscsys.h>
31 #endif
32 
33 DEAL_II_NAMESPACE_OPEN
34 
35 
36 // forward declarations
37 template <int dim> class Quadrature;
38 
39 
40 template <typename number> class Vector;
41 template <typename number> class FullMatrix;
42 template <typename number> class SparseMatrix;
43 
44 template <typename number> class BlockSparseMatrix;
45 template <typename Number> class BlockVector;
46 
47 template <int dim, int spacedim> class Mapping;
48 template <int dim, int spacedim> class DoFHandler;
49 
50 namespace hp
51 {
52  template <int> class QCollection;
53  template <int, int> class MappingCollection;
54  template <int, int> class DoFHandler;
55 }
56 
57 
58 #ifdef DEAL_II_WITH_PETSC
59 namespace PETScWrappers
60 {
61  class MatrixBase;
62  class VectorBase;
63  namespace MPI
64  {
65  class BlockSparseMatrix;
66  class BlockVector;
67  }
68 }
69 #endif
70 
71 #ifdef DEAL_II_WITH_TRILINOS
72 namespace TrilinosWrappers
73 {
74  class SparseMatrix;
75  class BlockSparseMatrix;
76  namespace MPI
77  {
78  class Vector;
79  class BlockVector;
80  }
81 }
82 #endif
83 
84 
217 namespace MatrixCreator
218 {
238  template <int dim, int spacedim, typename number>
239  void create_mass_matrix (const Mapping<dim, spacedim> &mapping,
240  const DoFHandler<dim,spacedim> &dof,
241  const Quadrature<dim> &q,
242  SparseMatrix<number> &matrix,
243  const Function<spacedim,number> *const a = nullptr,
244  const ConstraintMatrix &constraints = ConstraintMatrix());
245 
250  template <int dim, int spacedim, typename number>
252  const Quadrature<dim> &q,
253  SparseMatrix<number> &matrix,
254  const Function<spacedim,number> *const a = nullptr,
255  const ConstraintMatrix &constraints = ConstraintMatrix());
256 
276  template <int dim, int spacedim, typename number>
277  void create_mass_matrix (const Mapping<dim, spacedim> &mapping,
278  const DoFHandler<dim,spacedim> &dof,
279  const Quadrature<dim> &q,
281  const Function<spacedim,number> &rhs,
282  Vector<number> &rhs_vector,
283  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const a = nullptr,
284  const ConstraintMatrix &constraints = ConstraintMatrix());
285 
290  template <int dim, int spacedim, typename number>
292  const Quadrature<dim> &q,
294  const Function<spacedim,number> &rhs,
295  Vector<number> &rhs_vector,
296  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const a = nullptr,
297  const ConstraintMatrix &constraints = ConstraintMatrix());
298 
302  template <int dim, int spacedim, typename number>
304  const hp::DoFHandler<dim,spacedim> &dof,
305  const hp::QCollection<dim> &q,
306  SparseMatrix<number> &matrix,
307  const Function<spacedim,number> *const a = nullptr,
308  const ConstraintMatrix &constraints = ConstraintMatrix());
309 
313  template <int dim, int spacedim, typename number>
315  const hp::QCollection<dim> &q,
316  SparseMatrix<number> &matrix,
317  const Function<spacedim,number> *const a = nullptr,
318  const ConstraintMatrix &constraints = ConstraintMatrix());
319 
323  template <int dim, int spacedim, typename number>
325  const hp::DoFHandler<dim,spacedim> &dof,
326  const hp::QCollection<dim> &q,
328  const Function<spacedim,number> &rhs,
329  Vector<number> &rhs_vector,
330  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const a = nullptr,
331  const ConstraintMatrix &constraints = ConstraintMatrix());
332 
336  template <int dim, int spacedim, typename number>
338  const hp::QCollection<dim> &q,
340  const Function<spacedim,number> &rhs,
341  Vector<number> &rhs_vector,
342  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const a = nullptr,
343  const ConstraintMatrix &constraints = ConstraintMatrix());
344 
345 
371  template <int dim, int spacedim, typename number>
373  const DoFHandler<dim,spacedim> &dof,
374  const Quadrature<dim-1> &q,
376  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
377  Vector<number> &rhs_vector,
378  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
379  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const weight = 0,
380  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
381 
382 
387  template <int dim, int spacedim, typename number>
389  const Quadrature<dim-1> &q,
391  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
392  Vector<number> &rhs_vector,
393  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
394  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const a = nullptr,
395  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
396 
400  template <int dim, int spacedim, typename number>
402  const hp::DoFHandler<dim,spacedim> &dof,
403  const hp::QCollection<dim-1> &q,
405  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
406  Vector<number> &rhs_vector,
407  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
408  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const a = nullptr,
409  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
410 
414  template <int dim, int spacedim, typename number>
416  const hp::QCollection<dim-1> &q,
418  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
419  Vector<number> &rhs_vector,
420  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
421  const Function<spacedim,typename numbers::NumberTraits<number>::real_type> *const a = nullptr,
422  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
423 
443  template <int dim, int spacedim>
444  void create_laplace_matrix (const Mapping<dim, spacedim> &mapping,
445  const DoFHandler<dim,spacedim> &dof,
446  const Quadrature<dim> &q,
447  SparseMatrix<double> &matrix,
448  const Function<spacedim> *const a = nullptr,
449  const ConstraintMatrix &constraints = ConstraintMatrix());
450 
455  template <int dim, int spacedim>
457  const Quadrature<dim> &q,
458  SparseMatrix<double> &matrix,
459  const Function<spacedim> *const a = nullptr,
460  const ConstraintMatrix &constraints = ConstraintMatrix());
461 
480  template <int dim, int spacedim>
481  void create_laplace_matrix (const Mapping<dim, spacedim> &mapping,
482  const DoFHandler<dim,spacedim> &dof,
483  const Quadrature<dim> &q,
484  SparseMatrix<double> &matrix,
485  const Function<spacedim> &rhs,
486  Vector<double> &rhs_vector,
487  const Function<spacedim> *const a = nullptr,
488  const ConstraintMatrix &constraints = ConstraintMatrix());
489 
494  template <int dim, int spacedim>
496  const Quadrature<dim> &q,
497  SparseMatrix<double> &matrix,
498  const Function<spacedim> &rhs,
499  Vector<double> &rhs_vector,
500  const Function<spacedim> *const a = nullptr,
501  const ConstraintMatrix &constraints = ConstraintMatrix());
502 
507  template <int dim, int spacedim>
509  const hp::DoFHandler<dim,spacedim> &dof,
510  const hp::QCollection<dim> &q,
511  SparseMatrix<double> &matrix,
512  const Function<spacedim> *const a = nullptr,
513  const ConstraintMatrix &constraints = ConstraintMatrix());
514 
519  template <int dim, int spacedim>
521  const hp::QCollection<dim> &q,
522  SparseMatrix<double> &matrix,
523  const Function<spacedim> *const a = nullptr,
524  const ConstraintMatrix &constraints = ConstraintMatrix());
525 
530  template <int dim, int spacedim>
532  const hp::DoFHandler<dim,spacedim> &dof,
533  const hp::QCollection<dim> &q,
534  SparseMatrix<double> &matrix,
535  const Function<spacedim> &rhs,
536  Vector<double> &rhs_vector,
537  const Function<spacedim> *const a = nullptr,
538  const ConstraintMatrix &constraints = ConstraintMatrix());
539 
544  template <int dim, int spacedim>
546  const hp::QCollection<dim> &q,
547  SparseMatrix<double> &matrix,
548  const Function<spacedim> &rhs,
549  Vector<double> &rhs_vector,
550  const Function<spacedim> *const a = nullptr,
551  const ConstraintMatrix &constraints = ConstraintMatrix());
552 
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.");
565 }
566 
567 
568 
760 namespace MatrixTools
761 {
767  using namespace MatrixCreator;
768 
773  template <typename number>
774  void
775  apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
776  SparseMatrix<number> &matrix,
777  Vector<number> &solution,
778  Vector<number> &right_hand_side,
779  const bool eliminate_columns = true);
780 
786  template <typename number>
787  void
788  apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
790  BlockVector<number> &solution,
791  BlockVector<number> &right_hand_side,
792  const bool eliminate_columns = true);
793 
794 #ifdef DEAL_II_WITH_PETSC
795 
830  void
832  (const std::map<types::global_dof_index,PetscScalar> &boundary_values,
834  PETScWrappers::VectorBase &solution,
835  PETScWrappers::VectorBase &right_hand_side,
836  const bool eliminate_columns = true);
837 
841  void
842  apply_boundary_values (const std::map<types::global_dof_index,PetscScalar> &boundary_values,
845  PETScWrappers::MPI::BlockVector &right_hand_side,
846  const bool eliminate_columns = true);
847 
848 #endif
849 
850 #ifdef DEAL_II_WITH_TRILINOS
851 
884  void
885  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
888  TrilinosWrappers::MPI::Vector &right_hand_side,
889  const bool eliminate_columns = true);
890 
895  void
896  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
899  TrilinosWrappers::MPI::BlockVector &right_hand_side,
900  const bool eliminate_columns = true);
901 #endif
902 
921  template <typename number>
922  void
923  local_apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
924  const std::vector<types::global_dof_index> &local_dof_indices,
925  FullMatrix<number> &local_matrix,
926  Vector<number> &local_rhs,
927  const bool eliminate_columns);
928 
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.");
937 }
938 
939 
940 
941 DEAL_II_NAMESPACE_CLOSE
942 
943 #endif
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)
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:79
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())
unsigned int boundary_id
Definition: types.h:110
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.
Definition: dof_tools.h:46
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:335
static ::ExceptionBase & ExcComponentMismatch()
Definition: hp.h:102
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 >())