Flow123d  release_2.2.0-914-gf1a3a4f
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SchurComplement Class Reference

#include <schur.hh>

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

Public Member Functions

 SchurComplement (IS ia, Distribution *ds)
 
 SchurComplement (SchurComplement &other)
 
void set_tolerances (double r_tol, double a_tol, unsigned int max_it) override
 
void set_from_input (const Input::Record in_rec) override
 
LinSysget_system () const
 
Distributionget_distribution () const
 
 ~SchurComplement ()
 
void form_rhs ()
 
void set_complement (LinSys_PETSC *ls)
 Set complement LinSys object. More...
 
Distributionmake_complement_distribution ()
 get distribution of complement object if complement is defined More...
 
double get_solution_precision () override
 get precision of solving More...
 
LinSys::SolveInfo solve () override
 
void resolve ()
 
double compute_residual () override
 
- Public Member Functions inherited from LinSys_PETSC
 LinSys_PETSC (const Distribution *rows_ds, const std::string &params="")
 
 LinSys_PETSC (LinSys_PETSC &other)
 
const Distributionget_ds ()
 
const Mat * get_matrix () override
 
const Vec * get_rhs () override
 
PetscErrorCode set_matrix (Mat &matrix, MatStructure str) override
 
PetscErrorCode set_rhs (Vec &rhs) override
 
PetscErrorCode mat_zero_entries () override
 
PetscErrorCode rhs_zero_entries () override
 
void start_allocation () override
 
void start_add_assembly () override
 
void start_insert_assembly () override
 
void mat_set_values (int nrow, int *rows, int ncol, int *cols, double *vals) override
 
void rhs_set_values (int nrow, int *rows, double *vals) override
 
void preallocate_values (int nrow, int *rows, int ncol, int *cols)
 
void preallocate_matrix ()
 
void finish_assembly () override
 
void finish_assembly (MatAssemblyType assembly_type)
 
void apply_constrains (double scalar=1.) override
 
void set_initial_guess_nonzero (bool set_nonzero=true)
 
double get_absolute_accuracy () override
 
void view () override
 
 ~LinSys_PETSC ()
 
- Public Member Functions inherited from LinSys
 LinSys (const Distribution *rows_ds)
 
 LinSys (LinSys &other)
 
unsigned int size ()
 
unsigned int vec_lsize ()
 
void set_matrix_changed ()
 
void set_rhs_changed ()
 
bool is_matrix_changed ()
 
bool is_rhs_changed ()
 
const Vec & get_solution ()
 
void set_solution (double *sol_array)
 
double * get_solution_array ()
 
void mat_set_value (int row, int col, double val)
 
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_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)
 
double get_residual_norm ()
 
double get_relative_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 ~LinSys ()
 

Protected Member Functions

void create_inversion_matrix ()
 create IA matrix More...
 
void form_schur ()
 

Protected Attributes

Mat IA
 
Mat B
 
Mat Bt
 
Mat C
 
Mat xA
 
Mat IAB
 
int loc_size_A
 
int loc_size_B
 
IS IsA
 
IS IsB
 
Vec RHS1
 
Vec RHS2
 
Vec Sol1
 
Vec Sol2
 
SchurState state
 
int orig_lsize
 Size of local vector part of original system. More...
 
LinSys_PETSCCompl
 
Distributionds_
 
- Protected Attributes inherited from LinSys_PETSC
std::string params_
 command-line-like options for the PETSc solver More...
 
bool init_guess_nonzero
 flag for starting from nonzero guess More...
 
Mat matrix_
 Petsc matrix of the problem. More...
 
Vec rhs_
 PETSc vector constructed with vx array. More...
 
Vec residual_
 
double * v_rhs_
 local RHS array pointing to Vec rhs_ More...
 
Vec on_vec_
 Vectors for counting non-zero entries in diagonal block. More...
 
Vec off_vec_
 Vectors for counting non-zero entries in off-diagonal block. More...
 
double solution_precision_
 
KSP system
 
KSPConvergedReason reason
 
- Protected Attributes inherited from LinSys
double r_tol_
 relative tolerance of linear solver More...
 
double a_tol_
 absolute tolerance 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_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_
 

Additional Inherited Members

- Public Types inherited from LinSys_PETSC
typedef LinSys FactoryBaseType
 
- Public Types inherited from LinSys
enum  SetValuesMode {
  INSERT =INSERT_VALUES, ADD =ADD_VALUES, ALLOCATE, DONE,
  NONE
}
 
- Static Public Member Functions inherited from LinSys_PETSC
static const Input::Type::Recordget_input_type ()
 
- Static Public Member Functions inherited from LinSys
static Input::Type::Abstractget_input_type ()
 
- Protected Types inherited from LinSys
typedef std::pair< unsigned, double > Constraint_
 
typedef std::vector< Constraint_ConstraintVec_
 
- Static Protected Attributes inherited from LinSys
static constexpr double default_r_tol_ = 1e-7
 
static constexpr double default_a_tol_ = 1e-11
 
static constexpr unsigned int default_max_it_ = 1000
 

Detailed Description

Definition at line 64 of file schur.hh.

Constructor & Destructor Documentation

SchurComplement::SchurComplement ( IS  ia,
Distribution ds 
)

Constructor

Gets linear system with original matrix A and creates its inversion (IA matrix)

In current implementation the index set IsA has to be continuous sequence at the beginning of the local block of indices.

Create Schur complement system.

Parameters
[in]orig: original system
[in]inv_a: inversion of the A block
[in]ia: index set of the A block, default continuous given by inv_a: proc 1 2 3

Orig: ****** ****** **** IA : *** ** ***

Definition at line 61 of file schur.cc.

SchurComplement::SchurComplement ( SchurComplement other)

Copy constructor.

Definition at line 89 of file schur.cc.

SchurComplement::~SchurComplement ( )

Destructor. In particular it also delete complement linear system if it was passed in through the set_complement() method.

SCHUR COMPLEMENT destructor

Definition at line 398 of file schur.cc.

Member Function Documentation

double SchurComplement::compute_residual ( )
overridevirtual

The solve or resolve must be called prior to computing the residual.

Reimplemented from LinSys_PETSC.

Definition at line 386 of file schur.cc.

void SchurComplement::create_inversion_matrix ( )
protected

create IA matrix

Definition at line 270 of file schur.cc.

Here is the caller graph for this function:

void SchurComplement::form_rhs ( )

Compute only right hand side. This is useful when you change only rhs of the original system. TODO: We should ask original system if the matrix has changed (using LazyDependency) and possibly call only form_rhs, then this can be protected

Definition at line 233 of file schur.cc.

Here is the caller graph for this function:

void SchurComplement::form_schur ( )
protected

COMPUTE A SCHUR COMPLEMENT OF A PETSC MATRIX

given symmetric original matrix Orig has form A B x_1 RHS_1 B' C * x_2 = RHS_2 where the first block is given by index set IsA, and the second block by IsB user has to provide inverse IA of the A-block we suppose that original matrix have non-zero pattern for the schur complement

we return: Shur - schur complement, ShurRHS - RHS of the complemented system: (B' * IA * B - C) * x_2 = (B' * IA * RHS_1 - RHS_2) IAB - a matrix to compute eliminated part of the solution: x_1 = IA * RHS_1 - IAB * x_2

Actually as B' is stored separetly, the routine can be used also for nonsymetric original system

Definition at line 138 of file schur.cc.

Here is the caller graph for this function:

Distribution* SchurComplement::get_distribution ( ) const
inline

Returns distribution of the original system (solved by class SchurComplement).

Definition at line 95 of file schur.hh.

double SchurComplement::get_solution_precision ( )
overridevirtual

get precision of solving

Reimplemented from LinSys_PETSC.

Definition at line 339 of file schur.cc.

LinSys* SchurComplement::get_system ( ) const
inline

Returns pointer to LinSys object representing the schur complement.

Definition at line 90 of file schur.hh.

Distribution * SchurComplement::make_complement_distribution ( )

get distribution of complement object if complement is defined

Definition at line 264 of file schur.cc.

void SchurComplement::resolve ( )

Only resolve the system with current solution vector. This is necessary for nonlinear solvers.

  • If matrix and/or rhs is changed the Schur complement is formed.
  • Resolve is called to reconstruct eliminated part of the solution vector.

COMPUTE ELIMINATED PART OF THE ORIG. SYS. & RESTORE RHS and SOLUTION VECTORS x_1 = IA * RHS_1 - IAB * x_2

Definition at line 374 of file schur.cc.

Here is the caller graph for this function:

void SchurComplement::set_complement ( LinSys_PETSC ls)

Set complement LinSys object.

Definition at line 257 of file schur.cc.

void SchurComplement::set_from_input ( const Input::Record  in_rec)
overridevirtual

Sets specific parameters defined by user in input file and used to calculate. Call set_from_input of complement

Reimplemented from LinSys_PETSC.

Definition at line 110 of file schur.cc.

void SchurComplement::set_tolerances ( double  r_tol,
double  a_tol,
unsigned int  max_it 
)
overridevirtual

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.

Reimplemented from LinSys_PETSC.

Definition at line 251 of file schur.cc.

LinSys::SolveInfo SchurComplement::solve ( )
overridevirtual

Solve the system.

  • If matrix and/or rhs is changed the Schur complement is formed.
  • The inner linear solver is called for the Schur complement
  • Resolve is called to reconstruct eliminated part of the solution vector.

Reimplemented from LinSys_PETSC.

Definition at line 348 of file schur.cc.

Member Data Documentation

Mat SchurComplement::B
protected

Definition at line 147 of file schur.hh.

Mat SchurComplement::Bt
protected

Definition at line 147 of file schur.hh.

Mat SchurComplement::C
protected

Definition at line 148 of file schur.hh.

LinSys_PETSC* SchurComplement::Compl
protected

Definition at line 159 of file schur.hh.

Distribution* SchurComplement::ds_
protected

Definition at line 161 of file schur.hh.

Mat SchurComplement::IA
protected

Definition at line 145 of file schur.hh.

Mat SchurComplement::IAB
protected

Definition at line 150 of file schur.hh.

IS SchurComplement::IsA
protected

Definition at line 152 of file schur.hh.

IS SchurComplement::IsB
protected

Definition at line 152 of file schur.hh.

int SchurComplement::loc_size_A
protected

Definition at line 151 of file schur.hh.

int SchurComplement::loc_size_B
protected

Definition at line 151 of file schur.hh.

int SchurComplement::orig_lsize
protected

Size of local vector part of original system.

Definition at line 157 of file schur.hh.

Vec SchurComplement::RHS1
protected

Definition at line 153 of file schur.hh.

Vec SchurComplement::RHS2
protected

Definition at line 153 of file schur.hh.

Vec SchurComplement::Sol1
protected

Definition at line 154 of file schur.hh.

Vec SchurComplement::Sol2
protected

Definition at line 154 of file schur.hh.

SchurState SchurComplement::state
protected

Definition at line 156 of file schur.hh.

Mat SchurComplement::xA
protected

Definition at line 149 of file schur.hh.


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