Flow123d  master-f44eb46
Classes | Public Types | Public Member Functions | Static Public Member Functions | Protected Types | Protected Attributes | Static Protected Attributes | Private Member Functions | Private Attributes | Friends | List of all members
LinSys Class Referenceabstract

#include <la_linsys_new.hh>

Inheritance diagram for LinSys:
Inheritance graph
[legend]
Collaboration diagram for LinSys:
Collaboration graph
[legend]

Classes

struct  SolveInfo
 

Public Types

enum  SetValuesMode {
  INSERT =INSERT_VALUES , ADD =ADD_VALUES , ALLOCATE , DONE ,
  NONE
}
 

Public Member Functions

 LinSys (const Distribution *rows_ds)
 
 LinSys (LinSys &other)
 
unsigned int size ()
 
unsigned int vec_lsize ()
 
virtual const Mat * get_matrix ()
 
virtual const Vec * get_rhs ()
 
void set_matrix_changed ()
 
void set_rhs_changed ()
 
virtual void set_tolerances (double r_tol, double a_tol, double d_tol, unsigned int max_it)=0
 
bool is_matrix_changed ()
 
bool is_rhs_changed ()
 
virtual PetscErrorCode set_matrix (Mat &, MatStructure)
 
virtual PetscErrorCode set_rhs (Vec &)
 
virtual PetscErrorCode mat_zero_entries ()
 
virtual PetscErrorCode rhs_zero_entries ()
 
const Vec & get_solution ()
 
void set_solution (Vec sol_vec)
 
void set_solution (double *sol_array)
 
void set_solution ()
 
double * get_solution_array ()
 
virtual void start_allocation ()
 
virtual void start_add_assembly ()
 
virtual void start_insert_assembly ()
 
virtual void finish_assembly ()=0
 
virtual void mat_set_values (int nrow, int *rows, int ncol, int *cols, double *vals)=0
 
void mat_set_value (int row, int col, double val)
 
virtual void rhs_set_values (int nrow, int *rows, double *vals)=0
 
void rhs_set_value (int row, double val)
 
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. More...
 
void set_local_system (LocalSystem &local)
 
void set_local_system (LocalSystem &local, const std::vector< LongIdx > &local_to_global_map)
 
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)
 
void add_constraint (int row, double value)
 
virtual void apply_constrains (double scalar)=0
 
virtual SolveInfo solve ()=0
 
double get_residual_norm ()
 
virtual double compute_residual ()=0
 
double get_relative_accuracy ()
 
virtual double get_absolute_accuracy ()
 
void set_symmetric (bool flag=true)
 
bool is_symmetric ()
 
void set_positive_definite (bool flag=true)
 
void set_negative_definite (bool flag=true)
 
bool is_positive_definite ()
 
bool is_negative_definite ()
 
bool is_new ()
 
bool is_preallocated ()
 
void set_spd_via_symmetric_general (bool flag=true)
 
bool is_spd_via_symmetric_general ()
 
virtual void view (string)
 
virtual void set_from_input (const Input::Record in_rec)
 
virtual double get_solution_precision ()=0
 
virtual ~LinSys ()
 

Static Public Member Functions

static Input::Type::Abstractget_input_type ()
 

Protected Types

typedef std::pair< unsigned, double > Constraint_
 
typedef std::vector< Constraint_ConstraintVec_
 

Protected Attributes

double r_tol_
 relative tolerance of linear solver More...
 
double a_tol_
 absolute tolerance of linear solver More...
 
double d_tol_
 tolerance for divergence of linear solver More...
 
unsigned int max_it_
 maximum number of iterations of linear solver More...
 
MPI_Comm comm_
 
SetValuesMode status_
 Set value status of the linear system. More...
 
const unsigned lsize_
 local number of matrix rows (non-overlapping division of rows) More...
 
unsigned size_
 global number of matrix rows, i.e. problem size More...
 
const Distributionrows_ds_
 final distribution of rows of MH matrix More...
 
bool symmetric_
 
bool positive_definite_
 
bool negative_definite_
 
bool spd_via_symmetric_general_
 
bool matrix_changed_
 true if the matrix was changed since the last solve More...
 
bool rhs_changed_
 true if the right hand side was changed since the last solve More...
 
Vec solution_
 PETSc vector constructed with vb array. More...
 
double * v_solution_
 local solution array pointing into Vec solution_ More...
 
bool own_vec_
 Indicates if the solution vector has been allocated by this class. More...
 
bool own_solution_
 Indicates if the solution array has been allocated by this class. More...
 
double residual_norm_
 local solution array pointing into Vec solution_ More...
 
ConstraintVec_ constraints_
 
std::vector< double > globalSolution_
 global solution in numbering for linear system More...
 
Input::Record in_rec_
 

Static Protected Attributes

static constexpr double default_r_tol_ = 1e-7
 
static constexpr double default_a_tol_ = 1e-11
 
static constexpr double default_d_tol_ = 10000
 
static constexpr unsigned int default_max_it_ = 1000
 

Private Member Functions

LinSysblock (i, j) MatrixArray *matrix() virtual solve(Vector solution
 

Private Attributes

LinSys Vector RHS
 

Friends

class SchurComplement
 

Detailed Description

LinSys - matrix and a particular way how to compute solution for given RHS, i.e. this class can perform action of the matrix inverse

This is abstract class for members of possible hierarchical tree of the whole system.

Definition at line 169 of file la_linsys_new.hh.

Member Typedef Documentation

◆ Constraint_

typedef std::pair<unsigned,double> LinSys::Constraint_
protected

Definition at line 113 of file linsys.hh.

◆ ConstraintVec_

Definition at line 114 of file linsys.hh.

Member Enumeration Documentation

◆ SetValuesMode

Enumerator
INSERT 
ADD 
ALLOCATE 
DONE 
NONE 

Definition at line 96 of file linsys.hh.

Constructor & Destructor Documentation

◆ LinSys() [1/2]

LinSys::LinSys ( const Distribution rows_ds)
inline

Constructor. Constructor of abstract class should not be called directly, but is used for initialization of member common to all particular solvers.

Parameters
comm- MPI communicator

TODO: Vector solution_ is now initialized to NULL, but it should be rather allocated in the constructor instead of the method set_solution().

Definition at line 127 of file linsys.hh.

◆ LinSys() [2/2]

LinSys::LinSys ( LinSys other)
inline

Copy constructor.

Definition at line 147 of file linsys.hh.

◆ ~LinSys()

virtual LinSys::~LinSys ( )
inlinevirtual

Definition at line 652 of file linsys.hh.

Member Function Documentation

◆ add_constraint()

void LinSys::add_constraint ( int  row,
double  value 
)
inline

Adds Dirichlet constrain.

Parameters
row- global number of row that should be eliminated.
value- solution value at the given row

Definition at line 514 of file linsys.hh.

◆ apply_constrains()

virtual void LinSys::apply_constrains ( double  scalar)
pure virtual

Apply constrains to assembled matrix. Constrains are given by pairs: global row index, value. i.e. typedef pair<unsigned int, double> Constrain;

What is th meaning of ( const double factor ) form Cambridge code?

Implemented in LinSys_PETSC, and LinSys_BDDC.

◆ block()

LinSys* LinSys::block ( ,
 
)
private

◆ compute_residual()

virtual double LinSys::compute_residual ( )
pure virtual

Explicitly compute residual and its norm for current solution.

Implemented in SchurComplement, LinSys_PETSC, LinSys_PERMON, and LinSys_BDDC.

Here is the caller graph for this function:

◆ finish_assembly()

virtual void LinSys::finish_assembly ( )
pure virtual

Finish assembly of the whole system. For PETSC this should call MatEndAssembly with MAT_FINAL_ASSEMBLY

Implemented in LinSys_PETSC, and LinSys_BDDC.

Here is the caller graph for this function:

◆ get_absolute_accuracy()

virtual double LinSys::get_absolute_accuracy ( )
inlinevirtual

Returns information on absolute solver accuracy

Reimplemented in LinSys_PETSC, and LinSys_PERMON.

Definition at line 554 of file linsys.hh.

◆ get_input_type()

it::Abstract & LinSys::get_input_type ( )
static

Definition at line 29 of file linsys.cc.

Here is the caller graph for this function:

◆ get_matrix()

virtual const Mat* LinSys::get_matrix ( )
inlinevirtual

Returns PETSC matrix (only for PETSC solvers)

If matrix is changed, method set_matrix_changed() must be called. Example: @CODE MatDiagonalSet(schur->get_matrix(), new_diagonal, ADD_VALUES); schur->set_matrix_changed(); @ENDCODE

Reimplemented in LinSys_PETSC.

Definition at line 187 of file linsys.hh.

Here is the caller graph for this function:

◆ get_relative_accuracy()

double LinSys::get_relative_accuracy ( )
inline

Returns information on relative solver accuracy

Definition at line 547 of file linsys.hh.

◆ get_residual_norm()

double LinSys::get_residual_norm ( )
inline

Returns norm of the initial right hand side

Definition at line 535 of file linsys.hh.

◆ get_rhs()

virtual const Vec* LinSys::get_rhs ( )
inlinevirtual

Returns RHS vector (only for PETSC solvers)

If vector is changed, method set_rhs_changed() must be called. Example: @CODE VecScale(schur->get_rhs(), -1.0); schur->set_rhs_changed(); @ENDCODE

Reimplemented in LinSys_PETSC.

Definition at line 203 of file linsys.hh.

Here is the caller graph for this function:

◆ get_solution()

const Vec& LinSys::get_solution ( )
inline

Returns PETSC vector with solution. Underlying array can be provided on construction.

Definition at line 282 of file linsys.hh.

Here is the caller graph for this function:

◆ get_solution_array()

double* LinSys::get_solution_array ( )
inline

Returns PETSC subarray with solution. Underlying array can be provided on construction.

Definition at line 325 of file linsys.hh.

Here is the caller graph for this function:

◆ get_solution_precision()

virtual double LinSys::get_solution_precision ( )
pure virtual

Get precision of solving

Implemented in SchurComplement, LinSys_PETSC, LinSys_PERMON, and LinSys_BDDC.

Here is the caller graph for this function:

◆ is_matrix_changed()

bool LinSys::is_matrix_changed ( )
inline

Returns true if the system matrix has changed since the last solve.

Definition at line 233 of file linsys.hh.

◆ is_negative_definite()

bool LinSys::is_negative_definite ( )
inline

Definition at line 600 of file linsys.hh.

Here is the caller graph for this function:

◆ is_new()

bool LinSys::is_new ( )
inline

TODO: In fact we want to know if the matrix is already preallocated However to do this we need explicit finalisation of preallocating cycle.

Definition at line 605 of file linsys.hh.

◆ is_positive_definite()

bool LinSys::is_positive_definite ( )
inline

Definition at line 597 of file linsys.hh.

Here is the caller graph for this function:

◆ is_preallocated()

bool LinSys::is_preallocated ( )
inline

Definition at line 609 of file linsys.hh.

◆ is_rhs_changed()

bool LinSys::is_rhs_changed ( )
inline

Returns true if the system RHS has changed since the last solve.

Definition at line 239 of file linsys.hh.

◆ is_spd_via_symmetric_general()

bool LinSys::is_spd_via_symmetric_general ( )
inline

Definition at line 625 of file linsys.hh.

◆ is_symmetric()

bool LinSys::is_symmetric ( )
inline

Definition at line 570 of file linsys.hh.

◆ mat_set_value()

void LinSys::mat_set_value ( int  row,
int  col,
double  val 
)
inline

Shortcut for assembling just one element into the matrix. Similarly we can provide method accepting armadillo matrices.

Definition at line 369 of file linsys.hh.

◆ mat_set_values()

virtual void LinSys::mat_set_values ( int  nrow,
int *  rows,
int  ncol,
int *  cols,
double *  vals 
)
pure virtual

Assembly full rectangular submatrix into the system matrix. Should be virtual, implemented differently in particular solvers.

Implemented in LinSys_PETSC, and LinSys_BDDC.

Here is the caller graph for this function:

◆ mat_zero_entries()

virtual PetscErrorCode LinSys::mat_zero_entries ( )
inlinevirtual

Clears entries of the matrix

Reimplemented in LinSys_PETSC, and LinSys_BDDC.

Definition at line 264 of file linsys.hh.

Here is the caller graph for this function:

◆ rhs_set_value()

void LinSys::rhs_set_value ( int  row,
double  val 
)
inline

Shorcut for assembling just one element into RHS vector.

Definition at line 381 of file linsys.hh.

Here is the caller graph for this function:

◆ rhs_set_values()

virtual void LinSys::rhs_set_values ( int  nrow,
int *  rows,
double *  vals 
)
pure virtual

Set values of the system right-hand side. Should be virtual, implemented differently in particular solvers.

Implemented in LinSys_PETSC, and LinSys_BDDC.

Here is the caller graph for this function:

◆ rhs_zero_entries()

virtual PetscErrorCode LinSys::rhs_zero_entries ( )
inlinevirtual

Clears entries of the right-hand side

Reimplemented in LinSys_PETSC, and LinSys_BDDC.

Definition at line 273 of file linsys.hh.

Here is the caller graph for this function:

◆ set_from_input()

virtual void LinSys::set_from_input ( const Input::Record  in_rec)
inlinevirtual

Sets basic parameters of LinSys defined by user in input file and used to calculate

Reimplemented in SchurComplement, LinSys_PETSC, and LinSys_BDDC.

Definition at line 641 of file linsys.hh.

Here is the caller graph for this function:

◆ set_local_system() [1/2]

void LinSys::set_local_system ( LocalSystem local)
inline

Definition at line 392 of file linsys.hh.

◆ set_local_system() [2/2]

void LinSys::set_local_system ( LocalSystem local,
const std::vector< LongIdx > &  local_to_global_map 
)
inline

Sets local system, considering that the dof indices are locallized on current processor.

Parameters
local_to_global_map- maps the local dof indices to global ones

Definition at line 412 of file linsys.hh.

◆ set_matrix()

virtual PetscErrorCode LinSys::set_matrix ( Mat &  ,
MatStructure   
)
inlinevirtual

Sets PETSC matrix (only for PETSC solvers)

Reimplemented in LinSys_PETSC.

Definition at line 246 of file linsys.hh.

◆ set_matrix_changed()

void LinSys::set_matrix_changed ( )
inline

Sets matrix changed flag.

Definition at line 212 of file linsys.hh.

Here is the caller graph for this function:

◆ set_negative_definite()

void LinSys::set_negative_definite ( bool  flag = true)
inline

Provides user knowledge about negative definiteness.

Definition at line 588 of file linsys.hh.

Here is the caller graph for this function:

◆ set_positive_definite()

void LinSys::set_positive_definite ( bool  flag = true)
inline

Provides user knowledge about positive definiteness.

Definition at line 576 of file linsys.hh.

Here is the caller graph for this function:

◆ set_rhs()

virtual PetscErrorCode LinSys::set_rhs ( Vec &  )
inlinevirtual

Sets RHS vector (only for PETSC solvers)

Reimplemented in LinSys_PETSC.

Definition at line 255 of file linsys.hh.

◆ set_rhs_changed()

void LinSys::set_rhs_changed ( )
inline

Sets rhs changed flag (only for PETSC solvers)

Definition at line 218 of file linsys.hh.

Here is the caller graph for this function:

◆ set_solution() [1/3]

void LinSys::set_solution ( )
inline

Create PETSc solution

Definition at line 314 of file linsys.hh.

Here is the caller graph for this function:

◆ set_solution() [2/3]

void LinSys::set_solution ( double *  sol_array)
inline

Create PETSc solution

Definition at line 303 of file linsys.hh.

◆ set_solution() [3/3]

void LinSys::set_solution ( Vec  sol_vec)
inline

Set PETSc solution

Definition at line 290 of file linsys.hh.

Here is the caller graph for this function:

◆ set_spd_via_symmetric_general()

void LinSys::set_spd_via_symmetric_general ( bool  flag = true)
inline

Provides user knowledge about positive definiteness via symmetric general approach. This is useful for solving Darcy flow by mixed hybrid method, where blocks on subdomains are saddle point but interface among subdomains is only at the block of Lagrange multipliers and is symmetric positive definite. Problem condensed to interface can thus be solved by PCG method, although original problem is saddle point.

Definition at line 619 of file linsys.hh.

◆ set_symmetric()

void LinSys::set_symmetric ( bool  flag = true)
inline

Provides user knowledge about symmetry.

Definition at line 561 of file linsys.hh.

Here is the caller graph for this function:

◆ set_tolerances()

virtual void LinSys::set_tolerances ( double  r_tol,
double  a_tol,
double  d_tol,
unsigned int  max_it 
)
pure virtual

Set relative tolerance, absolute tolerance, and maximum number of iterations of the linear solver.

For each of these three parameters we first look for the value at user input if not set we use the value provided to this method and finally the default values are set by the call of this method in the constructor.

Implemented in SchurComplement, LinSys_PETSC, and LinSys_BDDC.

Here is the caller graph for this function:

◆ set_values() [1/2]

void LinSys::set_values ( int  nrow,
int *  rows,
int  ncol,
int *  cols,
PetscScalar *  mat_vals,
PetscScalar *  rhs_vals 
)
inline

Set values in the system matrix and values in the right-hand side vector on corresponding rows.

Definition at line 386 of file linsys.hh.

Here is the caller graph for this function:

◆ set_values() [2/2]

void LinSys::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 
)
inline

Add given dense matrix to preallocation. Shortcut to assembly into matrix and RHS in one call, possibly apply Dirichlet boundary conditions. row_dofs - are global indices of rows of dense matrix and rows of dense vector @rhs in global system col_dofs - are global indices of columns of the matrix, and possibly

Application of Dirichlet conditions: 1) Rows with negative dofs are set to zero. 2) Cols with negative dofs are eliminated. 3) If there are entries on global diagonal. We determine value K either from diagonal of local matrix, or (if it is zero) from diagonal average.

Caveats:

  • can not set dirichlet condition on zero dof
  • Armadillo stores matrix in column first form (Fortran like) which makes it not well suited for passing local matrices.

Definition at line 452 of file linsys.hh.

◆ size()

unsigned int LinSys::size ( )
inline

Returns global system size.

Definition at line 163 of file linsys.hh.

◆ solve()

virtual SolveInfo LinSys::solve ( )
pure virtual

Solve the system and return convergence reason.

Implemented in SchurComplement, LinSys_PETSC, LinSys_PERMON, and LinSys_BDDC.

Here is the caller graph for this function:

◆ start_add_assembly()

virtual void LinSys::start_add_assembly ( )
inlinevirtual

Switch linear system into adding assembly. (the only one supported by triplets ??)

Reimplemented in LinSys_PETSC.

Definition at line 341 of file linsys.hh.

Here is the caller graph for this function:

◆ start_allocation()

virtual void LinSys::start_allocation ( )
inlinevirtual

Switch linear system into allocating assembly. (only for PETSC_MPIAIJ_preallocate_by_assembly)

Reimplemented in LinSys_PETSC.

Definition at line 333 of file linsys.hh.

Here is the caller graph for this function:

◆ start_insert_assembly()

virtual void LinSys::start_insert_assembly ( )
inlinevirtual

Switch linear system into insert assembly. (not currently used)

Reimplemented in LinSys_PETSC.

Definition at line 349 of file linsys.hh.

◆ vec_lsize()

unsigned int LinSys::vec_lsize ( )
inline

Returns local system size. (for partitioning of solution vectors) for PETSC_MPIAIJ it is also partitioning of the matrix

Definition at line 172 of file linsys.hh.

◆ view()

virtual void LinSys::view ( string  )
inlinevirtual

Output the system in the Matlab format possibly with given ordering. Rather we shoud provide output operator <<, since it is more flexible.

Reimplemented in LinSys_PETSC, and LinSys_PERMON.

Definition at line 633 of file linsys.hh.

Friends And Related Function Documentation

◆ SchurComplement

friend class SchurComplement
friend

Definition at line 91 of file linsys.hh.

Member Data Documentation

◆ a_tol_

double LinSys::a_tol_
protected

absolute tolerance of linear solver

Definition at line 666 of file linsys.hh.

◆ comm_

MPI_Comm LinSys::comm_
protected

Definition at line 670 of file linsys.hh.

◆ constraints_

ConstraintVec_ LinSys::constraints_
protected

Definition at line 693 of file linsys.hh.

◆ d_tol_

double LinSys::d_tol_
protected

tolerance for divergence of linear solver

Definition at line 667 of file linsys.hh.

◆ default_a_tol_

constexpr double LinSys::default_a_tol_ = 1e-11
staticconstexprprotected

Definition at line 661 of file linsys.hh.

◆ default_d_tol_

constexpr double LinSys::default_d_tol_ = 10000
staticconstexprprotected

Definition at line 662 of file linsys.hh.

◆ default_max_it_

constexpr unsigned int LinSys::default_max_it_ = 1000
staticconstexprprotected

Definition at line 663 of file linsys.hh.

◆ default_r_tol_

constexpr double LinSys::default_r_tol_ = 1e-7
staticconstexprprotected

Definition at line 660 of file linsys.hh.

◆ globalSolution_

std::vector<double> LinSys::globalSolution_
protected

global solution in numbering for linear system

Definition at line 695 of file linsys.hh.

◆ in_rec_

Input::Record LinSys::in_rec_
protected

Definition at line 697 of file linsys.hh.

◆ lsize_

const unsigned LinSys::lsize_
protected

local number of matrix rows (non-overlapping division of rows)

Definition at line 673 of file linsys.hh.

◆ matrix_changed_

bool LinSys::matrix_changed_
protected

true if the matrix was changed since the last solve

Definition at line 683 of file linsys.hh.

◆ max_it_

unsigned int LinSys::max_it_
protected

maximum number of iterations of linear solver

Definition at line 668 of file linsys.hh.

◆ negative_definite_

bool LinSys::negative_definite_
protected

Definition at line 680 of file linsys.hh.

◆ own_solution_

bool LinSys::own_solution_
protected

Indicates if the solution array has been allocated by this class.

Definition at line 689 of file linsys.hh.

◆ own_vec_

bool LinSys::own_vec_
protected

Indicates if the solution vector has been allocated by this class.

Definition at line 688 of file linsys.hh.

◆ positive_definite_

bool LinSys::positive_definite_
protected

Definition at line 679 of file linsys.hh.

◆ r_tol_

double LinSys::r_tol_
protected

relative tolerance of linear solver

Definition at line 665 of file linsys.hh.

◆ residual_norm_

double LinSys::residual_norm_
protected

local solution array pointing into Vec solution_

Definition at line 691 of file linsys.hh.

◆ RHS

LinSys Vector LinSys::RHS
private

Definition at line 173 of file la_linsys_new.hh.

◆ rhs_changed_

bool LinSys::rhs_changed_
protected

true if the right hand side was changed since the last solve

Definition at line 684 of file linsys.hh.

◆ rows_ds_

const Distribution* LinSys::rows_ds_
protected

final distribution of rows of MH matrix

Definition at line 676 of file linsys.hh.

◆ size_

unsigned LinSys::size_
protected

global number of matrix rows, i.e. problem size

Definition at line 674 of file linsys.hh.

◆ solution_

Vec LinSys::solution_
protected

PETSc vector constructed with vb array.

Definition at line 686 of file linsys.hh.

◆ spd_via_symmetric_general_

bool LinSys::spd_via_symmetric_general_
protected

Definition at line 681 of file linsys.hh.

◆ status_

SetValuesMode LinSys::status_
protected

Set value status of the linear system.

Definition at line 671 of file linsys.hh.

◆ symmetric_

bool LinSys::symmetric_
protected

Definition at line 678 of file linsys.hh.

◆ v_solution_

double* LinSys::v_solution_
protected

local solution array pointing into Vec solution_

Definition at line 687 of file linsys.hh.


The documentation for this class was generated from the following files: