Flow123d
master-3768d5dec
|
Go to the documentation of this file.
75 #include "petscmath.h"
86 namespace Input {
namespace Type {
class Abstract; } }
133 int lsizeInt =
static_cast<int>( rows_ds->
lsize() );
136 size_ =
static_cast<unsigned>( sizeInt );
155 ASSERT_PERMANENT(
false ).error(
"Using copy constructor of LinSys is not allowed!");
188 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function get_matrix is not implemented for linsys type\n.");
204 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function get_rhs is not implemented for linsys type\n.");
227 virtual void set_tolerances(
double r_tol,
double a_tol,
unsigned int max_it) = 0;
247 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function set_matrix is not implemented for linsys type \n.");
256 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function set_rhs is not implemented for linsys type \n.");
265 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function mat_zero_entries is not implemented for linsys type \n.");
274 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function vec_zero_entries is not implemented for linsys type \n.");
296 VecRestoreArray(
solution_, &out_array );
334 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function start_allocation is not implemented for linsys type \n.");
342 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function start_add_assembly is not implemented for linsys type \n.");
350 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function start_insert_assembly is not implemented for linsys type \n.");
385 inline void set_values(
int nrow,
int *
rows,
int ncol,
int *
cols,PetscScalar *mat_vals, PetscScalar *rhs_vals)
420 for (
uint i=0; i<m; i++)
421 row_dofs[i]= local_to_global_map[local.
row_dofs[i]];
422 for (
uint i=0; i<n; i++)
423 col_dofs[i]= local_to_global_map[local.
col_dofs[i]];
458 bool negative_row =
false;
459 bool negative_col =
false;
461 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
462 if (row_dofs[l_row] < 0) {
464 tmp.col(l_row).zeros();
468 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
469 if (col_dofs[l_col] < 0) {
470 tmp_rhs -= matrix.col(l_col) * col_solution[l_col];
471 tmp.row(l_col).zeros();
476 if (negative_row && negative_col) {
478 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
479 if (row_dofs[l_row] < 0)
480 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
481 if (col_dofs[l_col] < 0 && row_dofs[l_row] == col_dofs[l_col]) {
482 double new_diagonal = fabs(matrix.at(l_row,l_col));
483 if (new_diagonal == 0.0) {
484 if (matrix.is_square()) {
485 new_diagonal = arma::sum( abs(matrix.diag())) / matrix.n_rows;
487 new_diagonal = arma::accu( abs(matrix) ) / matrix.n_elem;
490 tmp.at(l_col, l_row) = new_diagonal;
491 tmp_rhs(l_row) = new_diagonal * row_solution[l_row];
497 for(
int &row : row_dofs) row=abs(row);
500 for(
int &col : col_dofs) col=abs(col);
503 mat_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])),
504 col_dofs.size(),
const_cast<int *
>(&(col_dofs[0])), tmp.memptr() );
505 rhs_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])), tmp_rhs.memptr() );
529 virtual SolveInfo
solve()=0;
634 ASSERT_PERMANENT(
false )(
typeid(*this).name()).error(
"Function view is not implemented for linsys type %s \n.");
double get_residual_norm()
const Vec & get_solution()
SolveInfo(int cr, int nits)
bool is_spd_via_symmetric_general()
double get_relative_accuracy()
unsigned int lsize(int proc) const
get local size
double residual_norm_
local solution array pointing into Vec solution_
const unsigned lsize_
local number of matrix rows (non-overlapping division of rows)
void set_values(int nrow, int *rows, int ncol, int *cols, PetscScalar *mat_vals, PetscScalar *rhs_vals)
Set values in the system matrix and values in the right-hand side vector on corresponding rows.
virtual PetscErrorCode set_rhs(Vec &)
void set_symmetric(bool flag=true)
Support classes for parallel programing.
static constexpr bool value
virtual PetscErrorCode rhs_zero_entries()
void chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
virtual PetscErrorCode mat_zero_entries()
void set_spd_via_symmetric_general(bool flag=true)
virtual SolveInfo solve()=0
void set_negative_definite(bool flag=true)
static constexpr double default_r_tol_
void set_positive_definite(bool flag=true)
void set_matrix_changed()
virtual void start_allocation()
void add_constraint(int row, double value)
virtual double compute_residual()=0
ArmaMat< double, N, M > mat
Vec solution_
PETSc vector constructed with vb array.
virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it)=0
unsigned size_
global number of matrix rows, i.e. problem size
double a_tol_
absolute tolerance of linear solver
void rhs_set_value(int row, double val)
arma::mat matrix
local system matrix
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
virtual const Mat * get_matrix()
virtual double get_absolute_accuracy()
bool rhs_changed_
true if the right hand side was changed since the last solve
virtual void apply_constrains(double scalar)=0
SetValuesMode status_
Set value status of the linear system.
void set_solution(double *sol_array)
double * v_solution_
local solution array pointing into Vec solution_
static constexpr double default_a_tol_
void mat_set_value(int row, int col, double val)
bool own_vec_
Indicates if the solution vector has been allocated by this class.
std::vector< double > globalSolution_
global solution in numbering for linear system
arma::vec rhs
local system RHS
virtual void start_add_assembly()
const Distribution * rows_ds_
final distribution of rows of MH matrix
virtual void rhs_set_values(int nrow, int *rows, double *vals)=0
std::pair< unsigned, double > Constraint_
virtual void set_from_input(const Input::Record in_rec)
double r_tol_
relative tolerance of linear solver
virtual PetscErrorCode set_matrix(Mat &, MatStructure)
void set_local_system(LocalSystem &local, const std::vector< LongIdx > &local_to_global_map)
virtual void finish_assembly()=0
void set_values(std::vector< int > &row_dofs, std::vector< int > &col_dofs, const arma::mat &matrix, const arma::vec &rhs, const arma::vec &row_solution, const arma::vec &col_solution)
ConstraintVec_ constraints_
virtual const Vec * get_rhs()
bool matrix_changed_
true if the matrix was changed since the last solve
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
virtual double get_solution_precision()=0
bool own_solution_
Indicates if the solution array has been allocated by this class.
std::vector< Constraint_ > ConstraintVec_
void set_solution(Vec sol_vec)
void set_local_system(LocalSystem &local)
static Input::Type::Abstract & get_input_type()
virtual void view(string)
bool is_positive_definite()
double * get_solution_array()
LinSys(const Distribution *rows_ds)
static constexpr unsigned int default_max_it_
unsigned int max_it_
maximum number of iterations of linear solver
virtual void start_insert_assembly()
bool is_negative_definite()
bool spd_via_symmetric_general_
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0