Reference documentation for deal.II version GIT relicensing-478-g3275795167 2024-04-24 07:10: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
relaxation_block.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2010 - 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_relaxation_block_h
16#define dealii_relaxation_block_h
17
18#include <deal.II/base/config.h>
19
22
25#include <deal.II/lac/vector.h>
26
27#include <vector>
28
30
52template <typename MatrixType,
53 typename InverseNumberType = typename MatrixType::value_type,
54 typename VectorType = Vector<double>>
55class RelaxationBlock : protected PreconditionBlockBase<InverseNumberType>
56{
57private:
61 using number = typename MatrixType::value_type;
62
66 using value_type = InverseNumberType;
67
68public:
73
81 {
82 public:
87 const double relaxation = 1.,
88 const bool invert_diagonal = true,
89 const bool same_diagonal = false,
92 const double threshold = 0.,
93 VectorType *temp_ghost_vector = nullptr);
94
100
105
113
123
128
138 double threshold = 0.;
139
149 unsigned int kernel_size = 0;
150
171 std::vector<std::vector<unsigned int>> order;
172
181 mutable VectorType *temp_ghost_vector;
182
186 std::size_t
188 };
189
199 void
200 initialize(const MatrixType &A, const AdditionalData &parameters);
201
207 void
209
225 void
227
228protected:
237 void
238 do_step(VectorType &dst,
239 const VectorType &prev,
240 const VectorType &src,
241 const bool backward) const;
242
249 SmartPointer<const MatrixType,
252
259
260private:
264 void
265 block_kernel(const size_type block_begin, const size_type block_end);
266};
267
268
282template <typename MatrixType,
283 typename InverseNumberType = typename MatrixType::value_type,
284 typename VectorType = Vector<double>>
286 : public virtual Subscriptor,
287 protected RelaxationBlock<MatrixType, InverseNumberType, VectorType>
288{
289public:
293 // RelaxationBlockJacobi();
294
298 using number = typename MatrixType::value_type;
299
304 AdditionalData;
305
309 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::initialize;
310
314 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::clear;
315
319 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::size;
323 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::inverse;
327 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::
332 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::inverse_svd;
336 using PreconditionBlockBase<InverseNumberType>::log_statistics;
340 void
341 step(VectorType &dst, const VectorType &rhs) const;
342
346 void
347 Tstep(VectorType &dst, const VectorType &rhs) const;
348
353 void
354 vmult(VectorType &dst, const VectorType &rhs) const;
355
360 void
361 Tvmult(VectorType &dst, const VectorType &rhs) const;
362};
363
364
378template <typename MatrixType,
379 typename InverseNumberType = typename MatrixType::value_type,
380 typename VectorType = Vector<double>>
382 : public virtual Subscriptor,
383 protected RelaxationBlock<MatrixType, InverseNumberType, VectorType>
384{
385public:
389 // RelaxationBlockSOR();
390
394 using number = typename MatrixType::value_type;
395
400 AdditionalData;
401
405 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::initialize;
406
410 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::clear;
411
415 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::size;
419 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::inverse;
423 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::
428 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::inverse_svd;
432 using PreconditionBlockBase<InverseNumberType>::log_statistics;
436 void
437 step(VectorType &dst, const VectorType &rhs) const;
438
442 void
443 Tstep(VectorType &dst, const VectorType &rhs) const;
444
449 void
450 vmult(VectorType &dst, const VectorType &rhs) const;
451
456 void
457 Tvmult(VectorType &dst, const VectorType &rhs) const;
458};
459
460
474template <typename MatrixType,
475 typename InverseNumberType = typename MatrixType::value_type,
476 typename VectorType = Vector<double>>
478 : public virtual Subscriptor,
479 protected RelaxationBlock<MatrixType, InverseNumberType, VectorType>
480{
481public:
485 using number = typename MatrixType::value_type;
486
491 AdditionalData;
492
496 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::initialize;
497
501 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::clear;
502
506 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::size;
510 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::inverse;
514 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::
519 using RelaxationBlock<MatrixType, InverseNumberType, VectorType>::inverse_svd;
523 using PreconditionBlockBase<InverseNumberType>::log_statistics;
527 void
528 step(VectorType &dst, const VectorType &rhs) const;
529
533 void
534 Tstep(VectorType &dst, const VectorType &rhs) const;
535
540 void
541 vmult(VectorType &dst, const VectorType &rhs) const;
542
547 void
548 Tvmult(VectorType &dst, const VectorType &rhs) const;
549};
550
551
553
554#endif
Householder< number > & inverse_householder(size_type i)
FullMatrix< number > & inverse(size_type i)
LAPACKFullMatrix< number > & inverse_svd(size_type i)
void Tstep(VectorType &dst, const VectorType &rhs) const
void Tvmult(VectorType &dst, const VectorType &rhs) const
void vmult(VectorType &dst, const VectorType &rhs) const
void step(VectorType &dst, const VectorType &rhs) const
typename MatrixType::value_type number
void Tvmult(VectorType &dst, const VectorType &rhs) const
void step(VectorType &dst, const VectorType &rhs) const
void vmult(VectorType &dst, const VectorType &rhs) const
typename MatrixType::value_type number
void Tstep(VectorType &dst, const VectorType &rhs) const
void Tstep(VectorType &dst, const VectorType &rhs) const
void Tvmult(VectorType &dst, const VectorType &rhs) const
typename MatrixType::value_type number
void step(VectorType &dst, const VectorType &rhs) const
void vmult(VectorType &dst, const VectorType &rhs) const
AdditionalData(const double relaxation=1., const bool invert_diagonal=true, const bool same_diagonal=false, const typename PreconditionBlockBase< InverseNumberType >::Inversion inversion=PreconditionBlockBase< InverseNumberType >::gauss_jordan, const double threshold=0., VectorType *temp_ghost_vector=nullptr)
PreconditionBlockBase< InverseNumberType >::Inversion inversion
std::size_t memory_consumption() const
std::vector< std::vector< unsigned int > > order
void initialize(const MatrixType &A, const AdditionalData &parameters)
void block_kernel(const size_type block_begin, const size_type block_end)
void do_step(VectorType &dst, const VectorType &prev, const VectorType &src, const bool backward) const
SmartPointer< const AdditionalData, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > additional_data
SmartPointer< const MatrixType, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > A
typename MatrixType::value_type number
void invert_diagblocks()
InverseNumberType value_type
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
unsigned int global_dof_index
Definition types.h:81