ug4
|
operator implementation of the local Schur complement More...
#include <feti.h>
Classes | |
struct | StepConv |
Public Types | |
typedef TAlgebra | algebra_type |
typedef TAlgebra::matrix_type | matrix_type |
typedef TAlgebra::vector_type | vector_type |
![]() | |
typedef TAlgebra::vector_type | codomain_function_type |
Range space. | |
typedef TAlgebra::vector_type | domain_function_type |
Domain space. | |
![]() | |
typedef Y | codomain_function_type |
Range space. | |
typedef X | domain_function_type |
Domain space. | |
![]() | |
typedef TAlgebra | algebra_type |
type of algebra | |
typedef TAlgebra::matrix_type | matrix_type |
type of matrix | |
typedef TAlgebra::vector_type | vector_type |
type of vector | |
![]() | |
typedef TAlgebra::vector_type | vector_type |
type of vector | |
Public Member Functions | |
virtual void | apply (vector_type &f, const vector_type &u) |
virtual void | apply_sub (vector_type &f, const vector_type &u) |
solves the system | |
void | clear_total_itercnt_of_inner_solvers () |
int | get_total_itercnt_of_inner_solvers () |
virtual void | init () |
initializes the solver for operator A | |
void | init (const vector_type &u) |
implementation of the operator for the solution dependent initialization. | |
LocalSchurComplement () | |
constructor | |
virtual const char * | name () const |
name of solver | |
void | print_statistic_of_inner_solver (bool bPrintOnlyAverages) |
prints some convergence statistic of inner solvers | |
void | set_dirichlet_solver (SmartPtr< ILinearOperatorInverse< vector_type > > dirichletSolver) |
sets a sequential Dirichlet solver | |
void | set_feti_layouts (FetiLayouts< algebra_type > &fetiLayouts) |
sets the primal layouts | |
void | set_matrix (SmartPtr< MatrixOperator< matrix_type, vector_type > > A) |
set original matrix from which the local Schur complement is constructed | |
void | set_statistic_type (std::string type) |
sets statistic slot where next iterate should be counted | |
virtual | ~LocalSchurComplement () |
![]() | |
virtual void | prepare (TAlgebra::vector_type &u) |
default implementation for IOperator interface | |
virtual | ~ILinearOperator () |
virtual destructor | |
![]() | |
virtual void | apply (Y &d, const X &u)=0 |
computes the nonlinear mapping d := N(u) | |
virtual void | prepare (X &u)=0 |
prepares domain and codomain functions for application | |
virtual | ~IOperator () |
virtual destructor | |
![]() | |
SmartPtr< IDebugWriter< algebra_type > > | debug_writer () |
returns the debug writer | |
ConstSmartPtr< IDebugWriter< algebra_type > > | debug_writer () const |
bool | debug_writer_valid () const |
returns true if the debug writer is set | |
DebugWritingObject () | |
DebugWritingObject (const DebugWritingObject< algebra_type > &parent) | |
clone constructor | |
DebugWritingObject (SmartPtr< IDebugWriter< algebra_type > > spDebugWriter) | |
virtual void | set_debug (SmartPtr< IDebugWriter< algebra_type > > spDebugWriter) |
set debug writer | |
virtual | ~DebugWritingObject () |
virtual destructor | |
![]() | |
virtual void | set_debug (SmartPtr< IVectorDebugWriter< vector_type > > spDebugWriter) |
set debug writer | |
SmartPtr< IVectorDebugWriter< vector_type > > | vector_debug_writer () |
returns the debug writer | |
ConstSmartPtr< IVectorDebugWriter< vector_type > > | vector_debug_writer () const |
bool | vector_debug_writer_valid () const |
returns true if the debug writer is set | |
VectorDebugWritingObject () | |
VectorDebugWritingObject (SmartPtr< IVectorDebugWriter< vector_type > > spDebugWriter) | |
void | write_debug (const vector_type &vec, const char *filename) |
writing debug output for a vector (if debug writer set) | |
virtual | ~VectorDebugWritingObject () |
virtual destructor | |
Protected Attributes | |
int | m_applyCnt |
std::map< std::string, std::vector< StepConv > > | m_mvStepConv |
matrix_type * | m_pDirichletMatrix |
FetiLayouts< algebra_type > * | m_pFetiLayouts |
matrix_type * | m_pMatrix |
SmartPtr< MatrixOperator< matrix_type, vector_type > > | m_spDirichletOperator |
SmartPtr< ILinearOperatorInverse< vector_type > > | m_spDirichletSolver |
SmartPtr< MatrixOperator< matrix_type, vector_type > > | m_spOperator |
std::string | m_statType |
int | m_totalIterCntOfInnerSolvers |
![]() | |
SmartPtr< IDebugWriter< algebra_type > > | m_spDebugWriter |
Debug Writer. | |
![]() | |
SmartPtr< IVectorDebugWriter< vector_type > > | m_spVectorDebugWriter |
Debug Writer. | |
Additional Inherited Members | |
![]() | |
void | enter_debug_writer_section (const char *secDir) |
enters a debugging section | |
void | enter_debug_writer_section (std::string secDir) |
enters a debugging section | |
void | leave_debug_writer_section () |
leaves a debugging section | |
void | write_debug (const matrix_type &mat, const char *filename) |
write debug output for a matrix (if debug writer set) | |
void | write_debug (const matrix_type &mat, std::string name) |
write debug output for a matrix (if debug writer set) | |
![]() | |
void | enter_vector_debug_writer_section (const char *secDir) |
enters a debugging section | |
void | enter_vector_debug_writer_section (std::string secDir) |
enters a debugging section | |
void | leave_vector_debug_writer_section () |
leaves a debugging section | |
void | print_debugger_message (const char *msg) |
prints a debugger message (listing all the sections) | |
void | print_debugger_message (std::string msg) |
prints a debugger message (listing all the sections) | |
virtual void | write_debug (const vector_type &vec, std::string name) |
writing debug output for a vector (if debug writer set) | |
operator implementation of the local Schur complement
This operator is the application of the local Schur complement \(S_{\Delta}^{i}\). The underlying matrix must have at least two layouts. The first layout, layout level 0, will be used to describe subdomain internal interfaces (i.e. "pure" processor interfaces), all other layouts are used to identify the boundary, layout level 1: \( \Delta \) (edges of subdomains), layout level 2: \( \Pi \) (vertices of subdomains, a.k.a. "cross points") - will be constructed here, and the Schur complement is build w.r.t. to these variables.
typedef TAlgebra ug::LocalSchurComplement< TAlgebra >::algebra_type |
typedef TAlgebra::matrix_type ug::LocalSchurComplement< TAlgebra >::matrix_type |
typedef TAlgebra::vector_type ug::LocalSchurComplement< TAlgebra >::vector_type |
ug::LocalSchurComplement< TAlgebra >::LocalSchurComplement | ( | ) |
constructor
References new().
|
inlinevirtual |
|
virtual |
applies the Schur complement built from matrix operator set via 'set_matrix()' to 'u' and returns the result 'f := S times u'
Implements ug::ILinearOperator< TAlgebra::vector_type, TAlgebra::vector_type >.
References FETI_PROFILE_BEGIN, ug::LocalSchurComplement< TAlgebra >::StepConv::lastDef3b, name, ug::LocalSchurComplement< TAlgebra >::StepConv::numIter3b, pcl::ProcRank(), ug::PST_ADDITIVE, ug::PST_CONSISTENT, SPNULL, UG_LOG_ALL_PROCS, and UG_THROW.
|
virtual |
solves the system
Implements ug::ILinearOperator< TAlgebra::vector_type, TAlgebra::vector_type >.
References FETI_PROFILE_BEGIN.
|
inline |
|
inline |
|
virtual |
initializes the solver for operator A
This method must be called, before the apply() method can be invoked. It has to be called each time, when the matrix has been replaced. A deep copy of the matrix is then constructed and in this copy the rows belonging to the \(\Delta\) and \(\Pi\) unknowns are set to identity rows. This matrix is used in the solution of the local dirichlet problem.
Implements ug::ILinearOperator< TAlgebra::vector_type, TAlgebra::vector_type >.
References UG_THROW.
|
inlinevirtual |
implementation of the operator for the solution dependent initialization.
Implements ug::ILinearOperator< TAlgebra::vector_type, TAlgebra::vector_type >.
References ug::LocalSchurComplement< TAlgebra >::init().
Referenced by ug::LocalSchurComplement< TAlgebra >::init().
|
inlinevirtual |
name of solver
void ug::LocalSchurComplement< TAlgebra >::print_statistic_of_inner_solver | ( | bool | bPrintOnlyAverages | ) |
prints some convergence statistic of inner solvers
References pcl::ProcessCommunicator::allreduce(), pcl::NumProcs(), PCL_DT_DOUBLE, PCL_DT_INT, PCL_RO_MAX, PCL_RO_MIN, PCL_RO_SUM, and UG_LOG.
Referenced by ug::FETISolver< TAlgebra >::print_statistic_of_inner_solver().
|
inline |
sets a sequential Dirichlet solver
This method sets the Dirichlet Solver that is used to invert the inner matrix \(A_{II}\)
References ug::LocalSchurComplement< TAlgebra >::m_spDirichletSolver.
|
inline |
sets the primal layouts
References ug::LocalSchurComplement< TAlgebra >::m_pFetiLayouts.
|
inline |
set original matrix from which the local Schur complement is constructed
Using this method, the original matrix A is set. Given the matrix in the form \( A = \sum\limits_{p=1}^{N} A^{(p)}\) in additive form, with
\begin{align*} A^{(p)} \begin{pmatrix} A_{II}^{(p)} & A_{I \Delta}^{(p)} & A_{I \Pi}^{(p)} \\ A_{\Delta I}^{(p)} & A_{\Delta \Delta}^{(p)} & A_{\Delta \Pi}^{(p)} \\ A_{\Pi I}^{(p)} & A_{\Pi \Delta}^{(p)} & A_{\Pi \Pi}^{(p)} \end{pmatrix} \end{align*}
the local Schur complement is the processwise application of the operator
\begin{align*} S_{\Delta \Delta} = A_{\Delta \Delta} - \begin{pmatrix} A_{I \Delta}^T & A_{\Pi \Delta}^T \end{pmatrix} \begin{pmatrix} A_{II} & A_{I \Pi} \\ A_{\Pi I} & A_{\Pi \Pi} \end{pmatrix}^{-1} \begin{pmatrix} A_{I \Delta} \\ A_{\Pi \Delta} \end{pmatrix} * \end{align*}
References ug::LocalSchurComplement< TAlgebra >::m_spOperator.
|
inline |
sets statistic slot where next iterate should be counted
References ug::LocalSchurComplement< TAlgebra >::m_statType.
|
protected |
|
protected |
|
protected |
|
protected |
Referenced by ug::LocalSchurComplement< TAlgebra >::set_feti_layouts().
|
protected |
|
protected |
|
protected |
Referenced by ug::LocalSchurComplement< TAlgebra >::set_dirichlet_solver().
|
protected |
Referenced by ug::LocalSchurComplement< TAlgebra >::set_matrix().
|
protected |
Referenced by ug::LocalSchurComplement< TAlgebra >::set_statistic_type().
|
protected |