Reference documentation for deal.II version GIT relicensing-249-g48dc7357c7 2024-03-29 12:30:02+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
matrix_tools.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2023 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_matrix_tools_h
16#define dealii_matrix_tools_h
17
18
19#include <deal.II/base/config.h>
20
23
25
27
28#include <map>
29
30#ifdef DEAL_II_WITH_PETSC
31# include <petscsys.h>
32#endif
33
35
36
37// forward declarations
38#ifndef DOXYGEN
39template <int dim>
40class Quadrature;
41
42
43template <typename number>
44class Vector;
45template <typename number>
46class FullMatrix;
47template <typename number>
48class SparseMatrix;
49
50template <typename number>
52template <typename Number>
53class BlockVector;
54
55template <int dim, int spacedim>
56class Mapping;
57template <int dim, int spacedim>
59class DoFHandler;
60
61namespace hp
62{
63 template <int>
64 class QCollection;
65 template <int, int>
66 class MappingCollection;
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
98
306namespace MatrixTools
307{
313 using namespace MatrixCreator;
314
319 template <typename number>
320 void
322 const std::map<types::global_dof_index, number> &boundary_values,
323 SparseMatrix<number> &matrix,
324 Vector<number> &solution,
325 Vector<number> &right_hand_side,
326 const bool eliminate_columns = true);
327
333 template <typename number>
334 void
336 const std::map<types::global_dof_index, number> &boundary_values,
338 BlockVector<number> &solution,
339 BlockVector<number> &right_hand_side,
340 const bool eliminate_columns = true);
341
342#ifdef DEAL_II_WITH_PETSC
350 void
352 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
355 PETScWrappers::VectorBase &right_hand_side,
356 const bool eliminate_columns = true);
357
361 void
363 const std::map<types::global_dof_index, PetscScalar> &boundary_values,
366 PETScWrappers::MPI::BlockVector &right_hand_side,
367 const bool eliminate_columns = true);
368
369#endif
370
371#ifdef DEAL_II_WITH_TRILINOS
405 void
407 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
410 TrilinosWrappers::MPI::Vector &right_hand_side,
411 const bool eliminate_columns = true);
412
417 void
419 const std::map<types::global_dof_index, TrilinosScalar> &boundary_values,
422 TrilinosWrappers::MPI::BlockVector &right_hand_side,
423 const bool eliminate_columns = true);
424#endif
425
444 template <typename number>
445 void
447 const std::map<types::global_dof_index, number> &boundary_values,
448 const std::vector<types::global_dof_index> &local_dof_indices,
449 FullMatrix<number> &local_matrix,
450 Vector<number> &local_rhs,
451 const bool eliminate_columns);
452
457 "You are providing a matrix whose subdivision into "
458 "blocks in either row or column direction does not use "
459 "the same blocks sizes as the solution vector or "
460 "right hand side vectors, respectively.");
461} // namespace MatrixTools
462
463
464
466
467#endif
Abstract base class for mapping classes.
Definition mapping.h:316
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_CXX20_REQUIRES(condition)
Definition config.h:177
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
static ::ExceptionBase & ExcBlocksDontMatch()
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:494
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)
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:117