Flow123d  last_with_con_2.0.0-663-gd0e2296
linsys_PETSC.hh
Go to the documentation of this file.
1 /*!
2  *
3  * Copyright (C) 2015 Technical University of Liberec. All rights reserved.
4  *
5  * This program is free software; you can redistribute it and/or modify it under
6  * the terms of the GNU General Public License version 3 as published by the
7  * Free Software Foundation. (http://www.gnu.org/licenses/gpl-3.0.en.html)
8  *
9  * This program is distributed in the hope that it will be useful, but WITHOUT
10  * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
11  * FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.
12  *
13  *
14  * @file linsys_PETSC.hh
15  * @brief Solver based on the original PETSc solver using MPIAIJ matrix and succesive Schur complement construction
16  * @author Jakub Sistek
17  */
18 
19 #ifndef LA_LINSYS_PETSC_HH_
20 #define LA_LINSYS_PETSC_HH_
21 
22 // derived from base linsys
23 #include "la/linsys.hh"
24 
25 #include "la/distribution.hh"
28 
29 #include "petscksp.h"
30 
31 class LinSys_PETSC : public LinSys
32 {
33 
34 public:
36 
37  static const Input::Type::Record & get_input_type();
38 
39  LinSys_PETSC(const Distribution * rows_ds);
40 
41  /**
42  * Copy constructor.
43  */
44  LinSys_PETSC( LinSys_PETSC &other );
45 
46 
47  void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override;
48 
49  /**
50  * Returns whole Distribution class for distribution of the solution.
51  */
52  inline const Distribution* get_ds( )
53  {
54  return rows_ds_;
55  }
56 
57  const Mat *get_matrix() override
58  {
59  return &matrix_;
60  }
61 
62  const Vec *get_rhs() override
63  {
64  return &rhs_;
65  }
66 
67  PetscErrorCode set_matrix(Mat &matrix, MatStructure str) override
68  {
69  matrix_changed_ = true;
70  return MatCopy(matrix, matrix_, str);
71  }
72 
73  PetscErrorCode set_rhs(Vec &rhs) override
74  {
75  rhs_changed_ = true;
76  return VecCopy(rhs, rhs_);
77  }
78 
79  PetscErrorCode mat_zero_entries() override
80  {
81  matrix_changed_ = true;
82  return MatZeroEntries(matrix_);
83  }
84 
85  PetscErrorCode rhs_zero_entries() override
86  {
87  rhs_changed_ = true;
88  return VecSet(rhs_, 0);
89  }
90 
91  void start_allocation() override;
92 
93  void start_add_assembly() override;
94 
95  void start_insert_assembly() override;
96 
97  void mat_set_values( int nrow, int *rows, int ncol, int *cols, double *vals ) override;
98 
99  void rhs_set_values( int nrow, int *rows, double *vals ) override;
100 
101  void preallocate_values(int nrow,int *rows,int ncol,int *cols);
102 
103  void preallocate_matrix();
104 
105  void finish_assembly() override;
106 
107  void finish_assembly( MatAssemblyType assembly_type );
108 
109  void apply_constrains( double scalar = 1. ) override;
110 
111  void set_initial_guess_nonzero(bool set_nonzero = true);
112 
113  int solve() override;
114 
115  /**
116  * Returns information on absolute solver accuracy
117  */
118  inline double get_absolute_accuracy() override {
119  return a_tol_;
120  };
121 
122  void view( ) override;
123 
124  /**
125  * Sets specific parameters of LinSys_PETSC defined by user in input file and used to calculate
126  */
127  void set_from_input(const Input::Record in_rec) override;
128 
129  double get_solution_precision() override;
130 
131  double compute_residual() override;
132 
133  ~LinSys_PETSC( );
134 
135 private:
136  /// Registrar of class to factory
137  static const int registrar;
138 
139  // make a pointer to the data array out of a std::vector
140  template<typename T>
142  {
143  if ( array.size() ) return &(array[0]);
144  return PETSC_NULL;
145  }
146 
147  // PetscScalar to double casting functor
148  struct PetscScalar2Double_ : public std::unary_function< PetscScalar, double >
149  {
150  double operator()( PetscScalar arg )
151  {
152  return static_cast<double>( arg );
153  }
154  };
155 
156 protected:
157 
158  std::string params_; //!< command-line-like options for the PETSc solver
159 
160  bool init_guess_nonzero; //!< flag for starting from nonzero guess
161 
162  Mat matrix_; //!< Petsc matrix of the problem.
163  Vec rhs_; //!< PETSc vector constructed with vx array.
165 
166  double *v_rhs_; //!< local RHS array pointing to Vec rhs_
167 
168  Vec on_vec_; //!< Vectors for counting non-zero entries in diagonal block.
169  Vec off_vec_; //!< Vectors for counting non-zero entries in off-diagonal block.
170 
171 
172  double solution_precision_; // precision of KSP system solver
173 
174  KSP system;
175  KSPConvergedReason reason;
176 
177 
178 };
179 
180 #endif /* LA_LINSYS_PETSC_HH_ */
const Distribution * get_ds()
Definition: linsys_PETSC.hh:52
void apply_constrains(double scalar=1.) override
const Mat * get_matrix() override
Definition: linsys_PETSC.hh:57
bool matrix_changed_
true if the matrix was changed since the last solve
Definition: linsys.hh:595
double * v_rhs_
local RHS array pointing to Vec rhs_
void set_from_input(const Input::Record in_rec) override
double get_absolute_accuracy() override
Wrappers for linear systems based on MPIAIJ and MATIS format.
LinSys_PETSC(const Distribution *rows_ds)
Definition: linsys_PETSC.cc:56
int solve() override
Vec rhs_
PETSc vector constructed with vx array.
void start_allocation() override
double operator()(PetscScalar arg)
void preallocate_values(int nrow, int *rows, int ncol, int *cols)
void start_add_assembly() override
bool rhs_changed_
true if the right hand side was changed since the last solve
Definition: linsys.hh:596
double get_solution_precision() override
void start_insert_assembly() override
T * makePetscPointer_(std::vector< T > &array)
Accessor to the data with type Type::Record.
Definition: accessors.hh:292
const Distribution * rows_ds_
final distribution of rows of MH matrix
Definition: linsys.hh:588
double a_tol_
absolute tolerance of linear solver
Definition: linsys.hh:579
PetscErrorCode set_matrix(Mat &matrix, MatStructure str) override
Definition: linsys_PETSC.hh:67
PetscErrorCode mat_zero_entries() override
Definition: linsys_PETSC.hh:79
LinSys FactoryBaseType
Definition: linsys_PETSC.hh:35
Vec on_vec_
Vectors for counting non-zero entries in diagonal block.
void preallocate_matrix()
double compute_residual() override
PetscErrorCode set_rhs(Vec &rhs) override
Definition: linsys_PETSC.hh:73
std::string params_
command-line-like options for the PETSc solver
Support classes for parallel programing.
KSPConvergedReason reason
internal::NamedArg< char > arg(StringRef name, const T &arg)
Definition: format.h:3291
const Vec * get_rhs() override
Definition: linsys_PETSC.hh:62
Vec off_vec_
Vectors for counting non-zero entries in off-diagonal block.
void finish_assembly() override
void view() override
void rhs_set_values(int nrow, int *rows, double *vals) override
Abstract linear system class.
void set_initial_guess_nonzero(bool set_nonzero=true)
static const Input::Type::Record & get_input_type()
Definition: linsys_PETSC.cc:31
void set_tolerances(double r_tol, double a_tol, unsigned int max_it) override
Definition: linsys_PETSC.cc:85
Mat matrix_
Petsc matrix of the problem.
static const int registrar
Registrar of class to factory.
Record type proxy class.
Definition: type_record.hh:171
PetscErrorCode rhs_zero_entries() override
Definition: linsys_PETSC.hh:85
void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals) override
bool init_guess_nonzero
flag for starting from nonzero guess
double solution_precision_