109 int lsizeInt =
static_cast<int>( rows_ds->
lsize() );
112 size_ =
static_cast<unsigned>( sizeInt );
131 OLD_ASSERT(
false,
"Using copy constructor of LinSys is not allowed!");
164 OLD_ASSERT(
false,
"Function get_matrix is not implemented for linsys type %s \n.",
typeid(*this).name() );
180 OLD_ASSERT(
false,
"Function get_rhs is not implemented for linsys type %s \n.",
typeid(*this).name() );
203 virtual void set_tolerances(
double r_tol,
double a_tol,
unsigned int max_it) = 0;
221 virtual PetscErrorCode
set_matrix(Mat &matrix, MatStructure str)
223 OLD_ASSERT(
false,
"Function set_matrix is not implemented for linsys type %s \n.",
typeid(*this).name() );
232 OLD_ASSERT(
false,
"Function set_rhs is not implemented for linsys type %s \n.",
typeid(*this).name() );
241 OLD_ASSERT(
false,
"Function mat_zero_entries is not implemented for linsys type %s \n.",
typeid(*this).name() );
250 OLD_ASSERT(
false,
"Function vec_zero_entries is not implemented for linsys type %s \n.",
typeid(*this).name() );
266 if (sol_array == NULL) {
291 OLD_ASSERT(
false,
"Function start_allocation is not implemented for linsys type %s \n.",
typeid(*this).name() );
299 OLD_ASSERT(
false,
"Function start_add_assembly is not implemented for linsys type %s \n.",
typeid(*this).name() );
307 OLD_ASSERT(
false,
"Function start_insert_assembly is not implemented for linsys type %s \n.",
typeid(*this).name() );
319 virtual void mat_set_values(
int nrow,
int *rows,
int ncol,
int *cols,
double *vals)=0;
342 inline void set_values(
int nrow,
int *rows,
int ncol,
int *cols,PetscScalar *mat_vals, PetscScalar *rhs_vals)
366 const arma::mat &matrix,
const arma::vec &rhs,
367 const arma::vec &row_solution,
const arma::vec &col_solution)
370 arma::mat tmp = matrix.t();
371 arma::vec tmp_rhs = rhs;
372 bool negative_row =
false;
373 bool negative_col =
false;
375 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
376 if (row_dofs[l_row] < 0) {
378 tmp.col(l_row).zeros();
382 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
383 if (col_dofs[l_col] < 0) {
384 tmp_rhs -= matrix.col(l_col) * col_solution[l_col];
385 tmp.row(l_col).zeros();
390 if (negative_row && negative_col) {
392 for(
unsigned int l_row = 0; l_row < row_dofs.size(); l_row++)
393 if (row_dofs[l_row] < 0)
394 for(
unsigned int l_col = 0; l_col < col_dofs.size(); l_col++)
395 if (col_dofs[l_col] < 0 && row_dofs[l_row] == col_dofs[l_col]) {
396 double new_diagonal = fabs(matrix.at(l_row,l_col));
397 if (new_diagonal == 0.0) {
398 if (matrix.is_square()) {
399 new_diagonal = arma::sum( abs(matrix.diag())) / matrix.n_rows;
401 new_diagonal = arma::accu( abs(matrix) ) / matrix.n_elem;
404 tmp.at(l_col, l_row) = new_diagonal;
405 tmp_rhs(l_row) = new_diagonal * row_solution[l_row];
411 for(
int &row : row_dofs) row=abs(row);
414 for(
int &col : col_dofs) col=abs(col);
417 mat_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])),
418 col_dofs.size(),
const_cast<int *
>(&(col_dofs[0])), tmp.memptr() );
419 rhs_set_values(row_dofs.size(),
const_cast<int *
>(&(row_dofs[0])), tmp_rhs.memptr() );
443 virtual int solve()=0;
548 OLD_ASSERT(
false,
"Function view is not implemented for linsys type %s \n.",
typeid(*this).name() );
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()
virtual void rhs_set_values(int nrow, int *rows, double *vals)=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)
const unsigned lsize_
local number of matrix rows (non-overlapping division of rows)
static constexpr double default_r_tol_
ConstraintVec_ constraints_
void set_negative_definite(bool flag=true)
virtual void start_allocation()
virtual void finish_assembly()=0
virtual double get_solution_precision()=0
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
bool rhs_changed_
true if the right hand side was changed since the last solve
LinSys(const Distribution *rows_ds)
virtual void set_tolerances(double r_tol, double a_tol, unsigned int max_it)=0
virtual void apply_constrains(double scalar)=0
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()
static constexpr unsigned int default_max_it_
const Distribution * rows_ds_
final distribution of rows of MH matrix
double a_tol_
absolute tolerance of linear solver
static constexpr double default_a_tol_
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_
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
Abstract linear system class.
bool is_positive_definite()
void set_matrix_changed()
void set_positive_definite(bool flag=true)
virtual PetscErrorCode set_matrix(Mat &matrix, MatStructure str)
static Input::Type::Abstract & get_input_type()
virtual void mat_set_values(int nrow, int *rows, int ncol, int *cols, double *vals)=0
bool own_solution_
Indicates if the solution array has been allocated by this class.
virtual const Mat * get_matrix()
virtual double get_absolute_accuracy()
virtual double compute_residual()=0
Vec solution_
PETSc vector constructed with vb array.
bool is_spd_via_symmetric_general()
void rhs_set_value(int row, double val)
unsigned int lsize(int proc) const
get local size