Reference documentation for deal.II version 8.5.1
relaxation_block.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2010 - 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__relaxation_block_h
17 #define dealii__relaxation_block_h
18 
19 #include <deal.II/base/subscriptor.h>
20 #include <deal.II/base/smartpointer.h>
21 #include <deal.II/lac/vector.h>
22 #include <deal.II/lac/precondition_block_base.h>
23 #include <deal.II/lac/sparsity_pattern.h>
24 
25 #include <vector>
26 #include <set>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
53 template <typename MatrixType,
54  typename InverseNumberType = typename MatrixType::value_type,
55  typename VectorType = Vector<double> >
57  protected PreconditionBlockBase<InverseNumberType>
58 {
59 private:
63  typedef typename MatrixType::value_type number;
64 
68  typedef InverseNumberType value_type;
69 
70 public:
75 
82  class AdditionalData : public Subscriptor
83  {
84  public:
88  AdditionalData (const double relaxation = 1.,
89  const bool invert_diagonal = true,
90  const bool same_diagonal = false,
93  const double threshold = 0.,
94  VectorType *temp_ghost_vector = NULL);
95 
101 
105  double relaxation;
106 
114 
124 
129 
137  double threshold;
138 
159  std::vector<std::vector<unsigned int> > order;
160 
169  mutable VectorType *temp_ghost_vector;
170 
174  std::size_t memory_consumption() const;
175  };
176 
186  void initialize (const MatrixType &A,
187  const AdditionalData &parameters);
188 
194  void clear();
195 
199  bool empty () const;
200 
206  size_type j) const;
207 
223  void invert_diagblocks();
224 
225 protected:
234  void do_step (
235  VectorType &dst,
236  const VectorType &prev,
237  const VectorType &src,
238  const bool backward) const;
239 
247 
252 
253 private:
257  void block_kernel(const size_type block_begin, const size_type block_end);
258 };
259 
260 
276 template<typename MatrixType,
277  typename InverseNumberType = typename MatrixType::value_type,
278  typename VectorType = Vector<double> >
279 class RelaxationBlockJacobi : public virtual Subscriptor,
280  protected RelaxationBlock<MatrixType, InverseNumberType, VectorType>
281 {
282 public:
286 // RelaxationBlockJacobi();
287 
291  typedef typename MatrixType::value_type number;
292 
297 
302 
307 
335  void step (VectorType &dst, const VectorType &rhs) const;
336 
340  void Tstep (VectorType &dst, const VectorType &rhs) const;
341 
345  std::size_t memory_consumption() const;
346 };
347 
348 
364 template<typename MatrixType,
365  typename InverseNumberType = typename MatrixType::value_type,
366  typename VectorType = Vector<double> >
367 class RelaxationBlockSOR : public virtual Subscriptor,
368  protected RelaxationBlock<MatrixType, InverseNumberType, VectorType>
369 {
370 public:
374 // RelaxationBlockSOR();
375 
379  typedef typename MatrixType::value_type number;
380 
385 
390 
395 
423  void step (VectorType &dst, const VectorType &rhs) const;
424 
428  void Tstep (VectorType &dst, const VectorType &rhs) const;
429 };
430 
431 
447 template<typename MatrixType,
448  typename InverseNumberType = typename MatrixType::value_type,
449  typename VectorType = Vector<double> >
450 class RelaxationBlockSSOR : public virtual Subscriptor,
451  protected RelaxationBlock<MatrixType, InverseNumberType, VectorType>
452 {
453 public:
457  typedef typename MatrixType::value_type number;
458 
463 
468 
473 
478 
502  void step (VectorType &dst, const VectorType &rhs) const;
503 
507  void Tstep (VectorType &dst, const VectorType &rhs) const;
508 };
509 
510 
511 DEAL_II_NAMESPACE_CLOSE
512 
513 #endif
std::size_t memory_consumption() const
bool empty() const
value_type el(size_type i, size_type j) const
SmartPointer< const AdditionalData, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > additional_data
InverseNumberType value_type
PreconditionBlockBase< InverseNumberType >::Inversion inversion
MatrixType::value_type number
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=NULL)
void block_kernel(const size_type block_begin, const size_type block_end)
SmartPointer< const MatrixType, RelaxationBlock< MatrixType, InverseNumberType, VectorType > > A
void Tstep(VectorType &dst, const VectorType &rhs) const
unsigned int global_dof_index
Definition: types.h:88
void step(VectorType &dst, const VectorType &rhs) const
MatrixType::value_type number
MatrixType::value_type number
void Tstep(VectorType &dst, const VectorType &rhs) const
std::vector< std::vector< unsigned int > > order
void invert_diagblocks()
void Tstep(VectorType &dst, const VectorType &rhs) const
types::global_dof_index size_type
void step(VectorType &dst, const VectorType &rhs) const
MatrixType::value_type number
void initialize(const MatrixType &A, const AdditionalData &parameters)
void step(VectorType &dst, const VectorType &rhs) const
std::size_t memory_consumption() const
void do_step(VectorType &dst, const VectorType &prev, const VectorType &src, const bool backward) const