61 #include <boost/exception/detail/error_info_impl.hpp> 62 #include <boost/exception/info.hpp> 75 #include "petscmath.h" 86 namespace Input {
namespace Type {
class Abstract; } }
106 : converged_reason(cr), n_iterations(nits)
128 : comm_( rows_ds->get_comm() ), status_( NONE ), lsize_( rows_ds->lsize() ), rows_ds_(rows_ds),
129 symmetric_( false ), positive_definite_( false ), negative_definite_( false ),
130 spd_via_symmetric_general_( false ), solution_(NULL), v_solution_(NULL),
133 int lsizeInt =
static_cast<int>( rows_ds->
lsize() );
136 size_ =
static_cast<unsigned>( sizeInt );
138 r_tol_ = default_r_tol_;
139 a_tol_ = default_a_tol_;
140 max_it_ = default_max_it_;
147 : r_tol_(other.r_tol_), a_tol_(other.a_tol_), max_it_(other.max_it_), comm_(other.comm_), status_(other.status_),
148 lsize_( other.rows_ds_->lsize() ), size_(other.size_), rows_ds_(other.rows_ds_), symmetric_(other.symmetric_),
149 positive_definite_(other.positive_definite_), negative_definite_( other.negative_definite_ ),
150 spd_via_symmetric_general_(other.spd_via_symmetric_general_), matrix_changed_(other.matrix_changed_),
151 rhs_changed_(other.rhs_changed_), residual_norm_(other.residual_norm_), constraints_(other.constraints_),
152 globalSolution_(other.globalSolution_), in_rec_(other.in_rec_)
155 OLD_ASSERT(
false,
"Using copy constructor of LinSys is not allowed!");
188 OLD_ASSERT(
false,
"Function get_matrix is not implemented for linsys type %s \n.",
typeid(*this).name() );
204 OLD_ASSERT(
false,
"Function get_rhs is not implemented for linsys type %s \n.",
typeid(*this).name() );
212 { matrix_changed_ =
true;}
218 { rhs_changed_ =
true; }
227 virtual void set_tolerances(
double r_tol,
double a_tol,
unsigned int max_it) = 0;
233 {
return matrix_changed_;}
239 {
return rhs_changed_;}
245 virtual PetscErrorCode
set_matrix(Mat &matrix, MatStructure str)
247 OLD_ASSERT(
false,
"Function set_matrix is not implemented for linsys type %s \n.",
typeid(*this).name() );
256 OLD_ASSERT(
false,
"Function set_rhs is not implemented for linsys type %s \n.",
typeid(*this).name() );
265 OLD_ASSERT(
false,
"Function mat_zero_entries is not implemented for linsys type %s \n.",
typeid(*this).name() );
274 OLD_ASSERT(
false,
"Function vec_zero_entries is not implemented for linsys type %s \n.",
typeid(*this).name() );
290 if (sol_array == NULL) {
291 v_solution_ =
new double[ rows_ds_->lsize() + 1 ];
292 own_solution_ =
true;
295 v_solution_ = sol_array;
296 own_solution_ =
false;
299 ierr = VecCreateMPIWithArray( comm_,1, rows_ds_->lsize(), PETSC_DECIDE, v_solution_, &solution_ ); CHKERRV( ierr );
315 OLD_ASSERT(
false,
"Function start_allocation is not implemented for linsys type %s \n.",
typeid(*this).name() );
323 OLD_ASSERT(
false,
"Function start_add_assembly is not implemented for linsys type %s \n.",
typeid(*this).name() );
331 OLD_ASSERT(
false,
"Function start_insert_assembly is not implemented for linsys type %s \n.",
typeid(*this).name() );
337 virtual void finish_assembly( )=0;
343 virtual void mat_set_values(
int nrow,
int *rows,
int ncol,
int *cols,
double *vals)=0;
350 { mat_set_values(1,&row,1,&col,&val); }
356 virtual void rhs_set_values(
int nrow,
int *rows,
double *vals)=0;
362 { rhs_set_values(1,&row,&val); }
366 inline void set_values(
int nrow,
int *rows,
int ncol,
int *cols,PetscScalar *mat_vals, PetscScalar *rhs_vals)
368 mat_set_values(nrow, rows, ncol, cols, mat_vals);
369 rhs_set_values(nrow, rows, rhs_vals);
374 arma::mat tmp = local.
matrix.t();
382 mat_set_values(local.
matrix.n_rows, (
int *)(local.
row_dofs.memptr()),
386 rhs_set_values(local.
matrix.n_rows, (
int *)(local.
row_dofs.memptr()),
413 const arma::mat &matrix,
const arma::vec &rhs,
414 const arma::vec &row_solution,
const arma::vec &col_solution)
417 arma::mat tmp = matrix.t();
418 arma::vec tmp_rhs = rhs;
419 bool negative_row =
false;
420 bool negative_col =
false;
422 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
423 if (row_dofs[l_row] < 0) {
425 tmp.col(l_row).zeros();
429 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
430 if (col_dofs[l_col] < 0) {
431 tmp_rhs -= matrix.col(l_col) * col_solution[l_col];
432 tmp.row(l_col).zeros();
437 if (negative_row && negative_col) {
439 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
440 if (row_dofs[l_row] < 0)
441 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
442 if (col_dofs[l_col] < 0 && row_dofs[l_row] == col_dofs[l_col]) {
443 double new_diagonal = fabs(matrix.at(l_row,l_col));
444 if (new_diagonal == 0.0) {
445 if (matrix.is_square()) {
446 new_diagonal = arma::sum( abs(matrix.diag())) / matrix.n_rows;
448 new_diagonal = arma::accu( abs(matrix) ) / matrix.n_elem;
451 tmp.at(l_col, l_row) = new_diagonal;
452 tmp_rhs(l_row) = new_diagonal * row_solution[l_row];
458 for(
int &row : row_dofs) row=abs(row);
461 for(
int &col : col_dofs) col=abs(col);
464 mat_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])),
465 col_dofs.size(),
const_cast<int *
>(&(col_dofs[0])), tmp.memptr() );
466 rhs_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])), tmp_rhs.memptr() );
476 constraints_.push_back( Constraint_( static_cast<unsigned>( row ), value ) );
485 virtual void apply_constrains(
double scalar )=0;
496 return residual_norm_;
502 virtual double compute_residual() =0;
525 set_positive_definite(
false);
526 set_negative_definite(
false);
531 {
return symmetric_; }
538 positive_definite_ = flag;
541 set_negative_definite(
false);
550 negative_definite_ = flag;
553 set_positive_definite(
false);
558 {
return positive_definite_; }
561 {
return negative_definite_; }
566 return ( status_ == NONE );
570 return ( status_ == INSERT || status_ == ADD);
581 spd_via_symmetric_general_ = flag;
582 if (flag) set_symmetric();
586 {
return spd_via_symmetric_general_; }
595 OLD_ASSERT(
false,
"Function view is not implemented for linsys type %s \n.",
typeid(*this).name() );
604 set_tolerances(default_r_tol_, default_a_tol_, default_max_it_);
610 virtual double get_solution_precision() = 0;
615 if ( solution_ ) { ierr = VecDestroy(&solution_); CHKERRV( ierr ); }
616 if ( own_solution_ )
delete[] v_solution_;
621 static constexpr
double default_r_tol_ = 1e-7;
622 static constexpr
double default_a_tol_ = 1e-11;
623 static constexpr
unsigned int default_max_it_ = 1000;
std::vector< double > globalSolution_
global solution in numbering for linear system
SetValuesMode status_
Set value status of the linear system.
virtual void start_insert_assembly()
bool is_negative_definite()
bool spd_via_symmetric_general_
std::pair< unsigned, double > Constraint_
double get_residual_norm()
void set_symmetric(bool flag=true)
double * v_solution_
local solution array pointing into Vec solution_
virtual void set_from_input(const Input::Record in_rec)
bool matrix_changed_
true if the matrix was changed since the last solve
virtual void start_add_assembly()
virtual PetscErrorCode mat_zero_entries()
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)
const unsigned lsize_
local number of matrix rows (non-overlapping division of rows)
ConstraintVec_ constraints_
void set_negative_definite(bool flag=true)
virtual void start_allocation()
void add_constraint(int row, double value)
void set_spd_via_symmetric_general(bool flag=true)
unsigned size_
global number of matrix rows, i.e. problem size
arma::vec rhs
local system RHS
bool rhs_changed_
true if the right hand side was changed since the last solve
LinSys(const Distribution *rows_ds)
void set_local_system(LocalSystem &local)
static constexpr bool value
void mat_set_value(int row, int col, double val)
Global macros to enhance readability and debugging, general constants.
virtual PetscErrorCode set_rhs(Vec &rhs)
double get_relative_accuracy()
unsigned int max_it_
maximum number of iterations of linear solver
double * get_solution_array()
const Vec & get_solution()
const Distribution * rows_ds_
final distribution of rows of MH matrix
double a_tol_
absolute tolerance of linear solver
double residual_norm_
local solution array pointing into Vec solution_
void set_solution(double *sol_array)
virtual const Vec * get_rhs()
virtual PetscErrorCode rhs_zero_entries()
Support classes for parallel programing.
#define MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, comm)
std::vector< Constraint_ > ConstraintVec_
arma::mat matrix
local system matrix
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...
double r_tol_
relative tolerance of linear solver
bool is_positive_definite()
void set_matrix_changed()
static const Input::Type::Record & get_input_type()
void set_positive_definite(bool flag=true)
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
bool own_solution_
Indicates if the solution array has been allocated by this class.
virtual const Mat * get_matrix()
virtual double get_absolute_accuracy()
Vec solution_
PETSc vector constructed with vb array.
SolveInfo(int cr, int nits)
bool is_spd_via_symmetric_general()
void eliminate_solution()
void rhs_set_value(int row, double val)
unsigned int lsize(int proc) const
get local size