Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07: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
sparse_decomposition.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2002 - 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_sparse_decomposition_h
16#define dealii_sparse_decomposition_h
17
18#include <deal.II/base/config.h>
19
21
22#include <cmath>
23
25
105template <typename number>
106class SparseLUDecomposition : protected SparseMatrix<number>,
107 public virtual Subscriptor
108{
109protected:
117
118public:
123
128 virtual ~SparseLUDecomposition() override = 0;
129
134 virtual void
135 clear() override;
136
141 {
142 public:
146 explicit AdditionalData(const double strengthen_diagonal = 0.,
147 const unsigned int extra_off_diagonals = 0,
148 const bool use_previous_sparsity = false,
149 const SparsityPattern *use_this_sparsity = nullptr);
150
158
170
180
192 };
193
209 template <typename somenumber>
210 void
212 const AdditionalData parameters);
213
218 bool
219 empty() const;
220
227 m() const;
228
235 n() const;
236
243 template <class OutVector, class InVector>
244 void
245 vmult_add(OutVector &dst, const InVector &src) const;
246
254 template <class OutVector, class InVector>
255 void
256 Tvmult_add(OutVector &dst, const InVector &src) const;
257
262 virtual std::size_t
264
274 double,
275 << "The strengthening parameter " << arg1
276 << " is not greater or equal than zero!");
278protected:
283 template <typename somenumber>
284 void
286
293 virtual void
295
305 virtual number
306 get_strengthen_diagonal(const number rowsum, const size_type row) const;
307
312
318 std::vector<const size_type *> prebuilt_lower_bound;
319
323 void
325
326private:
337};
338
340//---------------------------------------------------------------------------
341
342#ifndef DOXYGEN
343
344template <typename number>
345inline number
347 const number /*rowsum*/,
348 const size_type /*row*/) const
349{
350 return strengthen_diagonal;
351}
352
353
354
355template <typename number>
356inline bool
358{
360}
361
362
363template <typename number>
366{
368}
369
370
371template <typename number>
374{
376}
377
378// Note: This function is required for full compatibility with
379// the LinearOperator class. ::MatrixInterfaceWithVmultAdd
380// picks up the vmult_add function in the protected SparseMatrix
381// base class.
382template <typename number>
383template <class OutVector, class InVector>
384inline void
386 const InVector &src) const
387{
388 OutVector tmp;
389 tmp.reinit(dst);
390 this->vmult(tmp, src);
391 dst += tmp;
392}
393
394// Note: This function is required for full compatibility with
395// the LinearOperator class. ::MatrixInterfaceWithVmultAdd
396// picks up the vmult_add function in the protected SparseMatrix
397// base class.
398template <typename number>
399template <class OutVector, class InVector>
400inline void
402 const InVector &src) const
403{
404 OutVector tmp;
405 tmp.reinit(dst);
406 this->Tvmult(tmp, src);
407 dst += tmp;
408}
409
410//---------------------------------------------------------------------------
411
412
413template <typename number>
415 const double strengthen_diag,
416 const unsigned int extra_off_diag,
417 const bool use_prev_sparsity,
418 const SparsityPattern *use_this_spars)
419 : strengthen_diagonal(strengthen_diag)
420 , extra_off_diagonals(extra_off_diag)
421 , use_previous_sparsity(use_prev_sparsity)
422 , use_this_sparsity(use_this_spars)
423{}
424
425
426#endif // DOXYGEN
427
429
430#endif // dealii_sparse_decomposition_h
AdditionalData(const double strengthen_diagonal=0., const unsigned int extra_off_diagonals=0, const bool use_previous_sparsity=false, const SparsityPattern *use_this_sparsity=nullptr)
virtual number get_strengthen_diagonal(const number rowsum, const size_type row) const
virtual void strengthen_diagonal_impl()
void Tvmult_add(OutVector &dst, const InVector &src) const
size_type n() const
size_type m() const
typename SparseMatrix< number >::size_type size_type
virtual ~SparseLUDecomposition() override=0
void vmult_add(OutVector &dst, const InVector &src) const
void copy_from(const SparseMatrix< somenumber > &matrix)
virtual std::size_t memory_consumption() const
std::vector< const size_type * > prebuilt_lower_bound
virtual void clear() override
void initialize(const SparseMatrix< somenumber > &matrix, const AdditionalData parameters)
SparsityPattern * own_sparsity
size_type n() const
size_type m() const
bool empty() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:502
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:503
static ::ExceptionBase & ExcInvalidStrengthening(double arg1)
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:516