Reference documentation for deal.II version 9.3.3
\(\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
24
26
27#include <map>
28
29#ifdef DEAL_II_WITH_PETSC
30# include <petscsys.h>
31#endif
32
34
35
36// forward declarations
37#ifndef DOXYGEN
38template <int dim>
39class Quadrature;
40
41
42template <typename number>
43class Vector;
44template <typename number>
45class FullMatrix;
46template <typename number>
47class SparseMatrix;
48
49template <typename number>
51template <typename Number>
52class BlockVector;
53
54template <int dim, int spacedim>
55class Mapping;
56template <int dim, int spacedim>
57class DoFHandler;
58
59namespace hp
60{
61 template <int>
62 class QCollection;
63 template <int, int>
64 class MappingCollection;
65 template <int, int>
66 class DoFHandler;
67} // namespace hp
68
69
70# ifdef DEAL_II_WITH_PETSC
71namespace PETScWrappers
72{
73 class MatrixBase;
74 class VectorBase;
75 namespace MPI
76 {
78 class BlockVector;
79 } // namespace MPI
80} // namespace PETScWrappers
81# endif
82
83# ifdef DEAL_II_WITH_TRILINOS
84namespace TrilinosWrappers
85{
86 class SparseMatrix;
88 namespace MPI
89 {
90 class Vector;
91 class BlockVector;
92 } // namespace MPI
93} // namespace TrilinosWrappers
94# endif
95#endif
96
97
231{
251 template <int dim, int spacedim, typename number>
252 void
254 const Mapping<dim, spacedim> & mapping,
255 const DoFHandler<dim, spacedim> & dof,
256 const Quadrature<dim> & q,
258 const Function<spacedim, number> *const a = nullptr,
260
265 template <int dim, int spacedim, typename number>
266 void
268 const DoFHandler<dim, spacedim> & dof,
269 const Quadrature<dim> & q,
271 const Function<spacedim, number> *const a = nullptr,
273
293 template <int dim, int spacedim, typename number>
294 void
296 const Mapping<dim, spacedim> & mapping,
297 const DoFHandler<dim, spacedim> & dof,
298 const Quadrature<dim> & q,
300 const Function<spacedim, number> & rhs,
301 Vector<number> & rhs_vector,
302 const Function<spacedim, number> *const a = nullptr,
304
309 template <int dim, int spacedim, typename number>
310 void
312 const DoFHandler<dim, spacedim> & dof,
313 const Quadrature<dim> & q,
315 const Function<spacedim, number> & rhs,
316 Vector<number> & rhs_vector,
317 const Function<spacedim, number> *const a = nullptr,
319
323 template <int dim, int spacedim, typename number>
324 void
327 const DoFHandler<dim, spacedim> & dof,
328 const hp::QCollection<dim> & q,
330 const Function<spacedim, number> *const a = nullptr,
332
336 template <int dim, int spacedim, typename number>
337 void
339 const DoFHandler<dim, spacedim> & dof,
340 const hp::QCollection<dim> & q,
342 const Function<spacedim, number> *const a = nullptr,
344
348 template <int dim, int spacedim, typename number>
349 void
352 const DoFHandler<dim, spacedim> & dof,
353 const hp::QCollection<dim> & q,
355 const Function<spacedim, number> & rhs,
356 Vector<number> & rhs_vector,
357 const Function<spacedim, number> *const a = nullptr,
359
363 template <int dim, int spacedim, typename number>
364 void
366 const DoFHandler<dim, spacedim> & dof,
367 const hp::QCollection<dim> & q,
369 const Function<spacedim, number> & rhs,
370 Vector<number> & rhs_vector,
371 const Function<spacedim, number> *const a = nullptr,
373
374
400 template <int dim, int spacedim, typename number>
401 void
403 const Mapping<dim, spacedim> & mapping,
404 const DoFHandler<dim, spacedim> &dof,
405 const Quadrature<dim - 1> & q,
407 const std::map<types::boundary_id, const Function<spacedim, number> *>
408 & boundary_functions,
409 Vector<number> & rhs_vector,
410 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
411 const Function<spacedim, number> *const weight = 0,
412 std::vector<unsigned int> component_mapping = {});
413
414
419 template <int dim, int spacedim, typename number>
420 void
422 const DoFHandler<dim, spacedim> &dof,
423 const Quadrature<dim - 1> & q,
425 const std::map<types::boundary_id, const Function<spacedim, number> *>
426 & boundary_functions,
427 Vector<number> & rhs_vector,
428 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
429 const Function<spacedim, number> *const a = nullptr,
430 std::vector<unsigned int> component_mapping = {});
431
435 template <int dim, int spacedim, typename number>
436 void
439 const DoFHandler<dim, spacedim> & dof,
440 const hp::QCollection<dim - 1> & q,
442 const std::map<types::boundary_id, const Function<spacedim, number> *>
443 & boundary_functions,
444 Vector<number> & rhs_vector,
445 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
446 const Function<spacedim, number> *const a = nullptr,
447 std::vector<unsigned int> component_mapping = {});
448
452 template <int dim, int spacedim, typename number>
453 void
455 const DoFHandler<dim, spacedim> &dof,
456 const hp::QCollection<dim - 1> & q,
458 const std::map<types::boundary_id, const Function<spacedim, number> *>
459 & boundary_functions,
460 Vector<number> & rhs_vector,
461 std::vector<types::global_dof_index> & dof_to_boundary_mapping,
462 const Function<spacedim, number> *const a = nullptr,
463 std::vector<unsigned int> component_mapping = {});
464
484 template <int dim, int spacedim>
485 void
487 const Mapping<dim, spacedim> & mapping,
488 const DoFHandler<dim, spacedim> &dof,
489 const Quadrature<dim> & q,
491 const Function<spacedim> *const a = nullptr,
493
498 template <int dim, int spacedim>
499 void
501 const DoFHandler<dim, spacedim> &dof,
502 const Quadrature<dim> & q,
504 const Function<spacedim> *const a = nullptr,
506
525 template <int dim, int spacedim>
526 void
528 const Mapping<dim, spacedim> & mapping,
529 const DoFHandler<dim, spacedim> &dof,
530 const Quadrature<dim> & q,
532 const Function<spacedim> & rhs,
533 Vector<double> & rhs_vector,
534 const Function<spacedim> *const a = nullptr,
536
541 template <int dim, int spacedim>
542 void
544 const DoFHandler<dim, spacedim> &dof,
545 const Quadrature<dim> & q,
547 const Function<spacedim> & rhs,
548 Vector<double> & rhs_vector,
549 const Function<spacedim> *const a = nullptr,
551
555 template <int dim, int spacedim>
556 void
559 const DoFHandler<dim, spacedim> & dof,
560 const hp::QCollection<dim> & q,
562 const Function<spacedim> *const a = nullptr,
564
568 template <int dim, int spacedim>
569 void
571 const DoFHandler<dim, spacedim> &dof,
572 const hp::QCollection<dim> & q,
574 const Function<spacedim> *const a = nullptr,
576
580 template <int dim, int spacedim>
581 void
584 const DoFHandler<dim, spacedim> & dof,
585 const hp::QCollection<dim> & q,
587 const Function<spacedim> & rhs,
588 Vector<double> & rhs_vector,
589 const Function<spacedim> *const a = nullptr,
591
595 template <int dim, int spacedim>
596 void
598 const DoFHandler<dim, spacedim> &dof,
599 const hp::QCollection<dim> & q,
601 const Function<spacedim> & rhs,
602 Vector<double> & rhs_vector,
603 const Function<spacedim> *const a = nullptr,
605
610 "You are providing either a right hand side function or a "
611 "coefficient with a number of vector components that is "
612 "inconsistent with the rest of the arguments. If you do "
613 "provide a coefficient or right hand side function, then "
614 "it either needs to have as many components as the finite "
615 "element in use, or only a single vector component. In "
616 "the latter case, the same value will be taken for "
617 "each vector component of the finite element.");
618} // namespace MatrixCreator
619
620
621
829namespace MatrixTools
830{
836 using namespace MatrixCreator;
837
842 template <typename number>
843 void
845 const std::map<types::global_dof_index, number> &boundary_values,
847 Vector<number> & solution,
848 Vector<number> & right_hand_side,
849 const bool eliminate_columns = true);
850
856 template <typename number>
857 void
859 const std::map<types::global_dof_index, number> &boundary_values,
861 BlockVector<number> & solution,
862 BlockVector<number> & right_hand_side,
863 const bool eliminate_columns = true);
864
865#ifdef DEAL_II_WITH_PETSC
901 void
903 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
905 PETScWrappers::VectorBase & solution,
906 PETScWrappers::VectorBase & right_hand_side,
907 const bool eliminate_columns = true);
908
912 void
914 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
917 PETScWrappers::MPI::BlockVector & right_hand_side,
918 const bool eliminate_columns = true);
919
920#endif
921
922#ifdef DEAL_II_WITH_TRILINOS
956 void
958 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
961 TrilinosWrappers::MPI::Vector & right_hand_side,
962 const bool eliminate_columns = true);
963
968 void
970 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
973 TrilinosWrappers::MPI::BlockVector & right_hand_side,
974 const bool eliminate_columns = true);
975#endif
976
995 template <typename number>
996 void
998 const std::map<types::global_dof_index, number> &boundary_values,
999 const std::vector<types::global_dof_index> & local_dof_indices,
1000 FullMatrix<number> & local_matrix,
1001 Vector<number> & local_rhs,
1002 const bool eliminate_columns);
1003
1008 "You are providing a matrix whose subdivision into "
1009 "blocks in either row or column direction does not use "
1010 "the same blocks sizes as the solution vector or "
1011 "right hand side vectors, respectively.");
1012} // namespace MatrixTools
1013
1014
1015
1017
1018#endif
Abstract base class for mapping classes.
Definition: mapping.h:304
Definition: vector.h:110
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
static ::ExceptionBase & ExcBlocksDontMatch()
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcComponentMismatch()
@ matrix
Contents is actually a 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={})
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 >())
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 >())
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
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: hp.h:118