61 #include <boost/exception/detail/error_info_impl.hpp> 62 #include <boost/exception/info.hpp> 75 #include "petscmath.h" 87 namespace Input {
namespace Type {
class Abstract; } }
107 : converged_reason(cr), n_iterations(nits)
129 : comm_( rows_ds->get_comm() ), status_( NONE ), lsize_( rows_ds->lsize() ), rows_ds_(rows_ds),
130 symmetric_( false ), positive_definite_( false ), negative_definite_( false ),
131 spd_via_symmetric_general_( false ), solution_(NULL), v_solution_(NULL),
134 int lsizeInt =
static_cast<int>( rows_ds->
lsize() );
137 size_ =
static_cast<unsigned>( sizeInt );
139 r_tol_ = default_r_tol_;
140 a_tol_ = default_a_tol_;
141 max_it_ = default_max_it_;
148 : r_tol_(other.r_tol_), a_tol_(other.a_tol_), max_it_(other.max_it_), comm_(other.comm_), status_(other.status_),
149 lsize_( other.rows_ds_->lsize() ), size_(other.size_), rows_ds_(other.rows_ds_), symmetric_(other.symmetric_),
150 positive_definite_(other.positive_definite_), negative_definite_( other.negative_definite_ ),
151 spd_via_symmetric_general_(other.spd_via_symmetric_general_), matrix_changed_(other.matrix_changed_),
152 rhs_changed_(other.rhs_changed_), residual_norm_(other.residual_norm_), constraints_(other.constraints_),
153 globalSolution_(other.globalSolution_), in_rec_(other.in_rec_)
156 OLD_ASSERT(
false,
"Using copy constructor of LinSys is not allowed!");
189 OLD_ASSERT(
false,
"Function get_matrix is not implemented for linsys type %s \n.",
typeid(*this).name() );
205 OLD_ASSERT(
false,
"Function get_rhs is not implemented for linsys type %s \n.",
typeid(*this).name() );
213 { matrix_changed_ =
true;}
219 { rhs_changed_ =
true; }
228 virtual void set_tolerances(
double r_tol,
double a_tol,
unsigned int max_it) = 0;
234 {
return matrix_changed_;}
240 {
return rhs_changed_;}
246 virtual PetscErrorCode
set_matrix(Mat &matrix, MatStructure str)
248 OLD_ASSERT(
false,
"Function set_matrix is not implemented for linsys type %s \n.",
typeid(*this).name() );
257 OLD_ASSERT(
false,
"Function set_rhs is not implemented for linsys type %s \n.",
typeid(*this).name() );
266 OLD_ASSERT(
false,
"Function mat_zero_entries is not implemented for linsys type %s \n.",
typeid(*this).name() );
275 OLD_ASSERT(
false,
"Function vec_zero_entries is not implemented for linsys type %s \n.",
typeid(*this).name() );
291 if (sol_array == NULL) {
292 v_solution_ =
new double[ rows_ds_->lsize() + 1 ];
293 own_solution_ =
true;
296 v_solution_ = sol_array;
297 own_solution_ =
false;
300 ierr = VecCreateMPIWithArray( comm_,1, rows_ds_->lsize(), PETSC_DECIDE, v_solution_, &solution_ ); CHKERRV( ierr );
316 OLD_ASSERT(
false,
"Function start_allocation is not implemented for linsys type %s \n.",
typeid(*this).name() );
324 OLD_ASSERT(
false,
"Function start_add_assembly is not implemented for linsys type %s \n.",
typeid(*this).name() );
332 OLD_ASSERT(
false,
"Function start_insert_assembly is not implemented for linsys type %s \n.",
typeid(*this).name() );
338 virtual void finish_assembly( )=0;
344 virtual void mat_set_values(
int nrow,
int *rows,
int ncol,
int *cols,
double *vals)=0;
351 { mat_set_values(1,&row,1,&col,&val); }
357 virtual void rhs_set_values(
int nrow,
int *rows,
double *vals)=0;
363 { rhs_set_values(1,&row,&val); }
367 inline void set_values(
int nrow,
int *rows,
int ncol,
int *cols,PetscScalar *mat_vals, PetscScalar *rhs_vals)
369 mat_set_values(nrow, rows, ncol, cols, mat_vals);
370 rhs_set_values(nrow, rows, rhs_vals);
375 arma::mat tmp = local.
matrix.t();
383 mat_set_values(local.
matrix.n_rows, (
int *)(local.
row_dofs.memptr()),
387 rhs_set_values(local.
matrix.n_rows, (
int *)(local.
row_dofs.memptr()),
414 const arma::mat &matrix,
const arma::vec &rhs,
415 const arma::vec &row_solution,
const arma::vec &col_solution)
418 arma::mat tmp = matrix.t();
419 arma::vec tmp_rhs = rhs;
420 bool negative_row =
false;
421 bool negative_col =
false;
423 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
424 if (row_dofs[l_row] < 0) {
426 tmp.col(l_row).zeros();
430 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
431 if (col_dofs[l_col] < 0) {
432 tmp_rhs -= matrix.col(l_col) * col_solution[l_col];
433 tmp.row(l_col).zeros();
438 if (negative_row && negative_col) {
440 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
441 if (row_dofs[l_row] < 0)
442 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
443 if (col_dofs[l_col] < 0 && row_dofs[l_row] == col_dofs[l_col]) {
444 double new_diagonal = fabs(matrix.at(l_row,l_col));
445 if (new_diagonal == 0.0) {
446 if (matrix.is_square()) {
447 new_diagonal = arma::sum( abs(matrix.diag())) / matrix.n_rows;
449 new_diagonal = arma::accu( abs(matrix) ) / matrix.n_elem;
452 tmp.at(l_col, l_row) = new_diagonal;
453 tmp_rhs(l_row) = new_diagonal * row_solution[l_row];
459 for(
int &row : row_dofs) row=abs(row);
462 for(
int &col : col_dofs) col=abs(col);
465 mat_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])),
466 col_dofs.size(),
const_cast<int *
>(&(col_dofs[0])), tmp.memptr() );
467 rhs_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])), tmp_rhs.memptr() );
477 constraints_.push_back( Constraint_( static_cast<unsigned>( row ), value ) );
486 virtual void apply_constrains(
double scalar )=0;
497 return residual_norm_;
503 virtual double compute_residual() =0;
526 set_positive_definite(
false);
527 set_negative_definite(
false);
532 {
return symmetric_; }
539 positive_definite_ = flag;
542 set_negative_definite(
false);
551 negative_definite_ = flag;
554 set_positive_definite(
false);
559 {
return positive_definite_; }
562 {
return negative_definite_; }
567 return ( status_ == NONE );
571 return ( status_ == INSERT || status_ == ADD);
582 spd_via_symmetric_general_ = flag;
583 if (flag) set_symmetric();
587 {
return spd_via_symmetric_general_; }
596 OLD_ASSERT(
false,
"Function view is not implemented for linsys type %s \n.",
typeid(*this).name() );
605 set_tolerances(default_r_tol_, default_a_tol_, default_max_it_);
611 virtual double get_solution_precision() = 0;
615 if ( solution_ ) {
chkerr(VecDestroy(&solution_)); }
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 chkerr(unsigned int ierr)
Replacement of new/delete operator in the spirit of xmalloc.
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