![]() |
|
#include <iterative_linear_solver.h>
Public Member Functions | |
DampedJacobi (const double &omega=2.0/3.0) | |
Empty constructor. More... | |
~DampedJacobi () | |
Empty destructor. More... | |
DampedJacobi (const DampedJacobi &)=delete | |
Broken copy constructor. More... | |
void | operator= (const DampedJacobi &)=delete |
Broken assignment operator. More... | |
void | clean_up_memory () |
Cleanup data that's stored for resolve (if any has been stored) More... | |
void | smoother_setup (DoubleMatrixBase *matrix_pt) |
Setup: Pass pointer to the matrix and store in cast form. More... | |
void | extract_diagonal_entries (DoubleMatrixBase *matrix_pt) |
Function to extract the diagonal entries from the matrix. More... | |
void | smoother_solve (const DoubleVector &rhs, DoubleVector &solution) |
void | solve (Problem *const &problem_pt, DoubleVector &result) |
void | solve (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &solution) |
void | resolve (const DoubleVector &rhs, DoubleVector &result) |
unsigned | iterations () const |
Number of iterations taken. More... | |
![]() | |
Smoother () | |
Empty constructor. More... | |
virtual | ~Smoother () |
Virtual empty destructor. More... | |
template<typename MATRIX > | |
void | check_validity_of_solve_helper_inputs (MATRIX *const &matrix_pt, const DoubleVector &rhs, DoubleVector &solution, const double &n_dof) |
![]() | |
IterativeLinearSolver () | |
IterativeLinearSolver (const IterativeLinearSolver &)=delete | |
Broken copy constructor. More... | |
void | operator= (const IterativeLinearSolver &)=delete |
Broken assignment operator. More... | |
virtual | ~IterativeLinearSolver () |
Destructor (empty) More... | |
Preconditioner *& | preconditioner_pt () |
Access function to preconditioner. More... | |
Preconditioner *const & | preconditioner_pt () const |
Access function to preconditioner (const version) More... | |
double & | tolerance () |
Access to convergence tolerance. More... | |
unsigned & | max_iter () |
Access to max. number of iterations. More... | |
void | enable_doc_convergence_history () |
Enable documentation of the convergence history. More... | |
void | disable_doc_convergence_history () |
Disable documentation of the convergence history. More... | |
void | open_convergence_history_file_stream (const std::string &file_name, const std::string &zone_title="") |
void | close_convergence_history_file_stream () |
Close convergence history output stream. More... | |
double | jacobian_setup_time () const |
double | linear_solver_solution_time () const |
return the time taken to solve the linear system More... | |
virtual double | preconditioner_setup_time () const |
returns the the time taken to setup the preconditioner More... | |
void | enable_setup_preconditioner_before_solve () |
Setup the preconditioner before the solve. More... | |
void | disable_setup_preconditioner_before_solve () |
Don't set up the preconditioner before the solve. More... | |
void | enable_error_after_max_iter () |
Throw an error if we don't converge within max_iter. More... | |
void | disable_error_after_max_iter () |
Don't throw an error if we don't converge within max_iter (default). More... | |
void | enable_iterative_solver_as_preconditioner () |
void | disable_iterative_solver_as_preconditioner () |
![]() | |
LinearSolver () | |
Empty constructor, initialise the member data. More... | |
LinearSolver (const LinearSolver &dummy)=delete | |
Broken copy constructor. More... | |
void | operator= (const LinearSolver &)=delete |
Broken assignment operator. More... | |
virtual | ~LinearSolver () |
Empty virtual destructor. More... | |
void | enable_doc_time () |
Enable documentation of solve times. More... | |
void | disable_doc_time () |
Disable documentation of solve times. More... | |
bool | is_doc_time_enabled () const |
Is documentation of solve times enabled? More... | |
bool | is_resolve_enabled () const |
Boolean flag indicating if resolves are enabled. More... | |
virtual void | enable_resolve () |
virtual void | disable_resolve () |
virtual void | solve (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result) |
virtual void | solve_transpose (Problem *const &problem_pt, DoubleVector &result) |
virtual void | solve_transpose (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &result) |
virtual void | solve_transpose (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result) |
virtual void | resolve_transpose (const DoubleVector &rhs, DoubleVector &result) |
virtual void | enable_computation_of_gradient () |
void | disable_computation_of_gradient () |
void | reset_gradient () |
void | get_gradient (DoubleVector &gradient) |
function to access the gradient, provided it has been computed More... | |
![]() | |
DistributableLinearAlgebraObject () | |
Default constructor - create a distribution. More... | |
DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix)=delete | |
Broken copy constructor. More... | |
void | operator= (const DistributableLinearAlgebraObject &)=delete |
Broken assignment operator. More... | |
virtual | ~DistributableLinearAlgebraObject () |
Destructor. More... | |
LinearAlgebraDistribution * | distribution_pt () const |
access to the LinearAlgebraDistribution More... | |
unsigned | nrow () const |
access function to the number of global rows. More... | |
unsigned | nrow_local () const |
access function for the num of local rows on this processor. More... | |
unsigned | nrow_local (const unsigned &p) const |
access function for the num of local rows on this processor. More... | |
unsigned | first_row () const |
access function for the first row on this processor More... | |
unsigned | first_row (const unsigned &p) const |
access function for the first row on this processor More... | |
bool | distributed () const |
distribution is serial or distributed More... | |
bool | distribution_built () const |
void | build_distribution (const LinearAlgebraDistribution *const dist_pt) |
void | build_distribution (const LinearAlgebraDistribution &dist) |
Private Member Functions | |
void | solve_helper (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &solution) |
Private Attributes | |
MATRIX * | Matrix_pt |
Pointer to the matrix. More... | |
Vector< double > | Matrix_diagonal |
Vector containing the diagonal entries of A. More... | |
bool | Resolving |
bool | Matrix_can_be_deleted |
unsigned | Iterations |
Number of iterations taken. More... | |
double | Omega |
Damping factor. More... | |
//////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// //////////////////////////////////////////////////////////////////// Damped Jacobi "solver" templated by matrix type. The "solver" exists in many different incarnations: It's an IterativeLinearSolver, and a Smoother, all of which use the same basic iteration.
|
inline |
Empty constructor.
References Eigen::internal::omega(), and oomph::DampedJacobi< MATRIX >::Omega.
|
inline |
Empty destructor.
References oomph::DampedJacobi< MATRIX >::clean_up_memory().
|
delete |
Broken copy constructor.
|
inlinevirtual |
Cleanup data that's stored for resolve (if any has been stored)
Reimplemented from oomph::LinearSolver.
References oomph::DampedJacobi< MATRIX >::Matrix_can_be_deleted, and oomph::DampedJacobi< MATRIX >::Matrix_pt.
Referenced by oomph::DampedJacobi< MATRIX >::~DampedJacobi().
|
inline |
Function to extract the diagonal entries from the matrix.
References i, oomph::DampedJacobi< MATRIX >::Matrix_diagonal, oomph::DampedJacobi< MATRIX >::Matrix_pt, OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.
Referenced by oomph::DampedJacobi< MATRIX >::smoother_setup(), and oomph::DampedJacobi< MATRIX >::solve().
|
inlinevirtual |
Number of iterations taken.
Implements oomph::IterativeLinearSolver.
References oomph::DampedJacobi< MATRIX >::Iterations.
|
delete |
Broken assignment operator.
|
inlinevirtual |
Re-solve the system defined by the last assembled Jacobian and the rhs vector specified here. Solution is returned in the vector result.
Reimplemented from oomph::LinearSolver.
References oomph::DampedJacobi< MATRIX >::Matrix_pt, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::DampedJacobi< MATRIX >::Resolving, and oomph::DampedJacobi< MATRIX >::solve().
|
inlinevirtual |
Setup: Pass pointer to the matrix and store in cast form.
Implements oomph::Smoother.
References oomph::DampedJacobi< MATRIX >::extract_diagonal_entries(), oomph::DampedJacobi< MATRIX >::Matrix_can_be_deleted, and oomph::DampedJacobi< MATRIX >::Matrix_pt.
|
inlinevirtual |
The smoother_solve function performs fixed number of iterations on the system A*result=rhs. The number of (smoothing) iterations is the same as the max. number of iterations in the underlying IterativeLinearSolver class.
Implements oomph::Smoother.
References oomph::DampedJacobi< MATRIX >::Matrix_pt, BiharmonicTestFunctions1::solution(), oomph::DampedJacobi< MATRIX >::solve_helper(), and oomph::Smoother::Use_as_smoother.
|
inlinevirtual |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.
Reimplemented from oomph::LinearSolver.
References oomph::DistributableLinearAlgebraObject::build_distribution(), oomph::DistributableLinearAlgebraObject::distribution_pt(), oomph::LinearSolver::Enable_resolve, oomph::DampedJacobi< MATRIX >::extract_diagonal_entries(), oomph::DampedJacobi< MATRIX >::Matrix_can_be_deleted, oomph::DampedJacobi< MATRIX >::Matrix_pt, oomph::DampedJacobi< MATRIX >::Resolving, BiharmonicTestFunctions1::solution(), oomph::DampedJacobi< MATRIX >::solve_helper(), and oomph::Smoother::Use_as_smoother.
|
virtual |
Use damped Jacobi iteration as an IterativeLinearSolver: This obtains the Jacobian matrix J and the residual vector r (needed for the Newton method) from the problem's get_jacobian function and returns the result of Jx=r.
/////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// /////////////////////////////////////////////////////////////////// Solver: Takes pointer to problem and returns the results vector which contains the solution of the linear system defined by the problem's fully assembled Jacobian and residual vector.
Implements oomph::LinearSolver.
References oomph::TerminateHelper::clean_up_memory(), oomph::Problem::communicator_pt(), f(), oomph::Problem::get_jacobian(), oomph::Problem::ndof(), oomph::oomph_info, and oomph::TimingHelpers::timer().
Referenced by oomph::DampedJacobi< MATRIX >::resolve().
|
private |
This is where the actual work is done – different implementations for different matrix types.
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.
References oomph::DoubleVector::build(), i, Global_Variables::Iterations, oomph::BlackBoxFDNewtonSolver::Max_iter, oomph::DoubleMatrixBase::multiply(), oomph::DoubleVector::norm(), oomph::DistributableLinearAlgebraObject::nrow(), oomph::SarahBL::Omega, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, oomph::oomph_info, oomph::DoubleMatrixBase::residual(), BiharmonicTestFunctions1::solution(), oomph::Global_string_for_annotation::string(), oomph::TimingHelpers::timer(), and plotDoE::x.
Referenced by oomph::DampedJacobi< MATRIX >::smoother_solve(), and oomph::DampedJacobi< MATRIX >::solve().
|
private |
Number of iterations taken.
Referenced by oomph::DampedJacobi< MATRIX >::iterations().
|
private |
Boolean flag to indicate if the matrix pointed to be Matrix_pt can be deleted.
Referenced by oomph::DampedJacobi< MATRIX >::clean_up_memory(), oomph::DampedJacobi< MATRIX >::smoother_setup(), and oomph::DampedJacobi< MATRIX >::solve().
|
private |
Vector containing the diagonal entries of A.
Referenced by oomph::DampedJacobi< MATRIX >::extract_diagonal_entries().
|
private |
Pointer to the matrix.
Referenced by oomph::DampedJacobi< MATRIX >::clean_up_memory(), oomph::DampedJacobi< MATRIX >::extract_diagonal_entries(), oomph::DampedJacobi< MATRIX >::resolve(), oomph::DampedJacobi< MATRIX >::smoother_setup(), oomph::DampedJacobi< MATRIX >::smoother_solve(), and oomph::DampedJacobi< MATRIX >::solve().
|
private |
Damping factor.
Referenced by oomph::DampedJacobi< MATRIX >::DampedJacobi().
|
private |
Boolean flag to indicate if the solve is done in re-solve mode, bypassing setup of matrix and preconditioner
Referenced by oomph::DampedJacobi< MATRIX >::resolve(), and oomph::DampedJacobi< MATRIX >::solve().