Reference documentation for deal.II version 8.5.1
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 #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 SparseMatrix;
62  class Vector;
63  namespace MPI
64  {
65  class SparseMatrix;
66  class BlockSparseMatrix;
67  class Vector;
68  class BlockVector;
69  }
70 }
71 #endif
72 
73 #ifdef DEAL_II_WITH_TRILINOS
74 namespace TrilinosWrappers
75 {
76  class SparseMatrix;
77  class Vector;
78  class BlockSparseMatrix;
79  class BlockVector;
80  namespace MPI
81  {
82  class Vector;
83  class BlockVector;
84  }
85 }
86 #endif
87 
88 
221 namespace MatrixCreator
222 {
239  template <int dim, int spacedim, typename number>
240  void create_mass_matrix (const Mapping<dim, spacedim> &mapping,
241  const DoFHandler<dim,spacedim> &dof,
242  const Quadrature<dim> &q,
243  SparseMatrix<number> &matrix,
244  const Function<spacedim,number> *const a = 0,
245  const ConstraintMatrix &constraints = ConstraintMatrix());
246 
251  template <int dim, int spacedim, typename number>
253  const Quadrature<dim> &q,
254  SparseMatrix<number> &matrix,
255  const Function<spacedim,number> *const a = 0,
256  const ConstraintMatrix &constraints = ConstraintMatrix());
257 
274  template <int dim, int spacedim, typename number>
275  void create_mass_matrix (const Mapping<dim, spacedim> &mapping,
276  const DoFHandler<dim,spacedim> &dof,
277  const Quadrature<dim> &q,
278  SparseMatrix<number> &matrix,
279  const Function<spacedim,number> &rhs,
280  Vector<number> &rhs_vector,
281  const Function<spacedim,number> *const a = 0,
282  const ConstraintMatrix &constraints = ConstraintMatrix());
283 
288  template <int dim, int spacedim, typename number>
290  const Quadrature<dim> &q,
291  SparseMatrix<number> &matrix,
292  const Function<spacedim,number> &rhs,
293  Vector<number> &rhs_vector,
294  const Function<spacedim,number> *const a = 0,
295  const ConstraintMatrix &constraints = ConstraintMatrix());
296 
300  template <int dim, int spacedim, typename number>
302  const hp::DoFHandler<dim,spacedim> &dof,
303  const hp::QCollection<dim> &q,
304  SparseMatrix<number> &matrix,
305  const Function<spacedim,number> *const a = 0,
306  const ConstraintMatrix &constraints = ConstraintMatrix());
307 
311  template <int dim, int spacedim, typename number>
313  const hp::QCollection<dim> &q,
314  SparseMatrix<number> &matrix,
315  const Function<spacedim,number> *const a = 0,
316  const ConstraintMatrix &constraints = ConstraintMatrix());
317 
321  template <int dim, int spacedim, typename number>
323  const hp::DoFHandler<dim,spacedim> &dof,
324  const hp::QCollection<dim> &q,
325  SparseMatrix<number> &matrix,
326  const Function<spacedim,number> &rhs,
327  Vector<number> &rhs_vector,
328  const Function<spacedim,number> *const a = 0,
329  const ConstraintMatrix &constraints = ConstraintMatrix());
330 
334  template <int dim, int spacedim, typename number>
336  const hp::QCollection<dim> &q,
337  SparseMatrix<number> &matrix,
338  const Function<spacedim,number> &rhs,
339  Vector<number> &rhs_vector,
340  const Function<spacedim,number> *const a = 0,
341  const ConstraintMatrix &constraints = ConstraintMatrix());
342 
343 
366  template <int dim, int spacedim, typename number>
368  const DoFHandler<dim,spacedim> &dof,
369  const Quadrature<dim-1> &q,
370  SparseMatrix<number> &matrix,
371  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
372  Vector<number> &rhs_vector,
373  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
374  const Function<spacedim,number> *const weight = 0,
375  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
376 
377 
382  template <int dim, int spacedim, typename number>
384  const Quadrature<dim-1> &q,
385  SparseMatrix<number> &matrix,
386  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
387  Vector<number> &rhs_vector,
388  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
389  const Function<spacedim,number> *const a = 0,
390  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
391 
395  template <int dim, int spacedim, typename number>
397  const hp::DoFHandler<dim,spacedim> &dof,
398  const hp::QCollection<dim-1> &q,
399  SparseMatrix<number> &matrix,
400  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
401  Vector<number> &rhs_vector,
402  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
403  const Function<spacedim,number> *const a = 0,
404  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
405 
409  template <int dim, int spacedim, typename number>
411  const hp::QCollection<dim-1> &q,
412  SparseMatrix<number> &matrix,
413  const std::map<types::boundary_id, const Function<spacedim,number>*> &boundary_functions,
414  Vector<number> &rhs_vector,
415  std::vector<types::global_dof_index> &dof_to_boundary_mapping,
416  const Function<spacedim,number> *const a = 0,
417  std::vector<unsigned int> component_mapping = std::vector<unsigned int>());
418 
435  template <int dim, int spacedim>
436  void create_laplace_matrix (const Mapping<dim, spacedim> &mapping,
437  const DoFHandler<dim,spacedim> &dof,
438  const Quadrature<dim> &q,
439  SparseMatrix<double> &matrix,
440  const Function<spacedim> *const a = 0,
441  const ConstraintMatrix &constraints = ConstraintMatrix());
442 
447  template <int dim, int spacedim>
449  const Quadrature<dim> &q,
450  SparseMatrix<double> &matrix,
451  const Function<spacedim> *const a = 0,
452  const ConstraintMatrix &constraints = ConstraintMatrix());
453 
469  template <int dim, int spacedim>
470  void create_laplace_matrix (const Mapping<dim, spacedim> &mapping,
471  const DoFHandler<dim,spacedim> &dof,
472  const Quadrature<dim> &q,
473  SparseMatrix<double> &matrix,
474  const Function<spacedim> &rhs,
475  Vector<double> &rhs_vector,
476  const Function<spacedim> *const a = 0,
477  const ConstraintMatrix &constraints = ConstraintMatrix());
478 
483  template <int dim, int spacedim>
485  const Quadrature<dim> &q,
486  SparseMatrix<double> &matrix,
487  const Function<spacedim> &rhs,
488  Vector<double> &rhs_vector,
489  const Function<spacedim> *const a = 0,
490  const ConstraintMatrix &constraints = ConstraintMatrix());
491 
496  template <int dim, int spacedim>
498  const hp::DoFHandler<dim,spacedim> &dof,
499  const hp::QCollection<dim> &q,
500  SparseMatrix<double> &matrix,
501  const Function<spacedim> *const a = 0,
502  const ConstraintMatrix &constraints = ConstraintMatrix());
503 
508  template <int dim, int spacedim>
510  const hp::QCollection<dim> &q,
511  SparseMatrix<double> &matrix,
512  const Function<spacedim> *const a = 0,
513  const ConstraintMatrix &constraints = ConstraintMatrix());
514 
519  template <int dim, int spacedim>
521  const hp::DoFHandler<dim,spacedim> &dof,
522  const hp::QCollection<dim> &q,
523  SparseMatrix<double> &matrix,
524  const Function<spacedim> &rhs,
525  Vector<double> &rhs_vector,
526  const Function<spacedim> *const a = 0,
527  const ConstraintMatrix &constraints = ConstraintMatrix());
528 
533  template <int dim, int spacedim>
535  const hp::QCollection<dim> &q,
536  SparseMatrix<double> &matrix,
537  const Function<spacedim> &rhs,
538  Vector<double> &rhs_vector,
539  const Function<spacedim> *const a = 0,
540  const ConstraintMatrix &constraints = ConstraintMatrix());
541 
546  "You are providing either a right hand side function or a "
547  "coefficient with a number of vector components that is "
548  "inconsistent with the rest of the arguments. If you do "
549  "provide a coefficient or right hand side function, then "
550  "it either needs to have as many components as the finite "
551  "element in use, or only a single vector component. In "
552  "the latter case, the same value will be taken for "
553  "each vector component of the finite element.");
554 }
555 
556 
557 
749 namespace MatrixTools
750 {
756  using namespace MatrixCreator;
757 
762  template <typename number>
763  void
764  apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
765  SparseMatrix<number> &matrix,
766  Vector<number> &solution,
767  Vector<number> &right_hand_side,
768  const bool eliminate_columns = true);
769 
775  template <typename number>
776  void
777  apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
779  BlockVector<number> &solution,
780  BlockVector<number> &right_hand_side,
781  const bool eliminate_columns = true);
782 
783 #ifdef DEAL_II_WITH_PETSC
784 
803  void
804  apply_boundary_values (const std::map<types::global_dof_index,PetscScalar> &boundary_values,
806  PETScWrappers::Vector &solution,
807  PETScWrappers::Vector &right_hand_side,
808  const bool eliminate_columns = true);
809 
829  void
830  apply_boundary_values (const std::map<types::global_dof_index,PetscScalar> &boundary_values,
832  PETScWrappers::MPI::Vector &solution,
833  PETScWrappers::MPI::Vector &right_hand_side,
834  const bool eliminate_columns = true);
835 
839  void
840  apply_boundary_values (const std::map<types::global_dof_index,PetscScalar> &boundary_values,
843  PETScWrappers::MPI::BlockVector &right_hand_side,
844  const bool eliminate_columns = true);
845 
846 #endif
847 
848 #ifdef DEAL_II_WITH_TRILINOS
849 
866  void
867  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
869  TrilinosWrappers::Vector &solution,
870  TrilinosWrappers::Vector &right_hand_side,
871  const bool eliminate_columns = true);
872 
877  void
878  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
881  TrilinosWrappers::BlockVector &right_hand_side,
882  const bool eliminate_columns = true);
883 
903  void
904  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
907  TrilinosWrappers::MPI::Vector &right_hand_side,
908  const bool eliminate_columns = true);
909 
914  void
915  apply_boundary_values (const std::map<types::global_dof_index,TrilinosScalar> &boundary_values,
918  TrilinosWrappers::MPI::BlockVector &right_hand_side,
919  const bool eliminate_columns = true);
920 #endif
921 
940  template <typename number>
941  void
942  local_apply_boundary_values (const std::map<types::global_dof_index,number> &boundary_values,
943  const std::vector<types::global_dof_index> &local_dof_indices,
944  FullMatrix<number> &local_matrix,
945  Vector<number> &local_rhs,
946  const bool eliminate_columns);
947 
952  "You are providing a matrix whose subdivision into "
953  "blocks in either row or column direction does not use "
954  "the same blocks sizes as the solution vector or "
955  "right hand side vectors, respectively.");
956 }
957 
958 
959 
960 DEAL_II_NAMESPACE_CLOSE
961 
962 #endif
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=std::vector< unsigned int >())
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:80
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=0, const ConstraintMatrix &constraints=ConstraintMatrix())
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.
Definition: dof_tools.h:46
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:553
static ::ExceptionBase & ExcComponentMismatch()
Definition: hp.h:102
static ::ExceptionBase & ExcBlocksDontMatch()
unsigned char boundary_id
Definition: types.h:110