Dense LU decomposition-based solve of linear system assembled via finite differencing of the residuals Vector. Even more inefficient than DenseLU but excellent sanity check! More...
#include <linear_solver.h>
Public Member Functions | |
FD_LU () | |
Constructor: empty. | |
FD_LU (const FD_LU &dummy)=delete | |
Broken copy constructor. | |
void | operator= (const FD_LU &)=delete |
Broken assignment operator. | |
void | solve (Problem *const &problem_pt, DoubleVector &result) |
Solver: Takes pointer to problem and returns the results Vector which contains the solution of the linear system defined by the problem's residual Vector (Jacobian computed by FD approx.) | |
void | solve (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &result) |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system. | |
void | solve (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result) |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system Call the broken base-class version. If you want this, please implement it. | |
![]() | |
DenseLU () | |
Constructor, initialise storage. | |
DenseLU (const DenseLU &dummy)=delete | |
Broken copy constructor. | |
void | operator= (const DenseLU &)=delete |
Broken assignment operator. | |
~DenseLU () | |
Destructor, clean up the stored LU factors. | |
double | jacobian_setup_time () const |
returns the time taken to assemble the jacobian matrix and residual vector | |
virtual double | linear_solver_solution_time () const |
return the time taken to solve the linear system (needs to be overloaded for each linear solver) | |
![]() | |
LinearSolver () | |
Empty constructor, initialise the member data. | |
LinearSolver (const LinearSolver &dummy)=delete | |
Broken copy constructor. | |
void | operator= (const LinearSolver &)=delete |
Broken assignment operator. | |
virtual | ~LinearSolver () |
Empty virtual destructor. | |
void | enable_doc_time () |
Enable documentation of solve times. | |
void | disable_doc_time () |
Disable documentation of solve times. | |
bool | is_doc_time_enabled () const |
Is documentation of solve times enabled? | |
bool | is_resolve_enabled () const |
Boolean flag indicating if resolves are enabled. | |
virtual void | enable_resolve () |
Enable resolve (i.e. store matrix and/or LU decomposition, say) Virtual so it can be overloaded to perform additional tasks. | |
virtual void | disable_resolve () |
Disable resolve (i.e. store matrix and/or LU decomposition, say) This function simply resets an internal flag. It's virtual so it can be overloaded to perform additional tasks such as cleaning up memory that is only required for the resolve. | |
virtual void | solve_transpose (Problem *const &problem_pt, DoubleVector &result) |
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 (broken virtual). | |
virtual void | solve_transpose (DoubleMatrixBase *const &matrix_pt, const DoubleVector &rhs, DoubleVector &result) |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system. | |
virtual void | solve_transpose (DoubleMatrixBase *const &matrix_pt, const Vector< double > &rhs, Vector< double > &result) |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system. | |
virtual void | resolve (const DoubleVector &rhs, DoubleVector &result) |
Resolve the system defined by the last assembled jacobian and the rhs vector. Solution is returned in the vector result. (broken virtual) | |
virtual void | resolve_transpose (const DoubleVector &rhs, DoubleVector &result) |
Solver: Resolve the system defined by the last assembled jacobian and the rhs vector. Solution is returned in the vector result. (broken virtual) | |
virtual void | enable_computation_of_gradient () |
function to enable the computation of the gradient required for the globally convergent Newton method | |
void | disable_computation_of_gradient () |
function to disable the computation of the gradient required for the globally convergent Newton method | |
void | reset_gradient () |
function to reset the size of the gradient before each Newton solve | |
void | get_gradient (DoubleVector &gradient) |
function to access the gradient, provided it has been computed | |
![]() | |
DistributableLinearAlgebraObject () | |
Default constructor - create a distribution. | |
DistributableLinearAlgebraObject (const DistributableLinearAlgebraObject &matrix)=delete | |
Broken copy constructor. | |
void | operator= (const DistributableLinearAlgebraObject &)=delete |
Broken assignment operator. | |
virtual | ~DistributableLinearAlgebraObject () |
Destructor. | |
LinearAlgebraDistribution * | distribution_pt () const |
access to the LinearAlgebraDistribution | |
unsigned | nrow () const |
access function to the number of global rows. | |
unsigned | nrow_local () const |
access function for the num of local rows on this processor. | |
unsigned | nrow_local (const unsigned &p) const |
access function for the num of local rows on this processor. | |
unsigned | first_row () const |
access function for the first row on this processor | |
unsigned | first_row (const unsigned &p) const |
access function for the first row on this processor | |
bool | distributed () const |
distribution is serial or distributed | |
bool | distribution_built () const |
if the communicator_pt is null then the distribution is not setup then false is returned, otherwise return true | |
void | build_distribution (const LinearAlgebraDistribution *const dist_pt) |
setup the distribution of this distributable linear algebra object | |
void | build_distribution (const LinearAlgebraDistribution &dist) |
setup the distribution of this distributable linear algebra object | |
Additional Inherited Members | |
![]() | |
void | factorise (DoubleMatrixBase *const &matrix_pt) |
Perform the LU decomposition of the matrix. | |
void | backsub (const DoubleVector &rhs, DoubleVector &result) |
Do the backsubstitution step to solve the system LU result = rhs. | |
void | backsub (const Vector< double > &rhs, Vector< double > &result) |
perform back substitution using Vector<double> | |
void | clean_up_memory () |
Clean up the stored LU factors. | |
![]() | |
void | clear_distribution () |
clear the distribution of this distributable linear algebra object | |
![]() | |
double | Jacobian_setup_time |
Jacobian setup time. | |
double | Solution_time |
Solution time. | |
int | Sign_of_determinant_of_matrix |
Sign of the determinant of the matrix (obtained during the LU decomposition) | |
![]() | |
bool | Enable_resolve |
Boolean that indicates whether the matrix (or its factors, in the case of direct solver) should be stored so that the resolve function can be used. | |
bool | Doc_time |
Boolean flag that indicates whether the time taken. | |
bool | Compute_gradient |
flag that indicates whether the gradient required for the globally convergent Newton method should be computed or not | |
bool | Gradient_has_been_computed |
flag that indicates whether the gradient was computed or not | |
DoubleVector | Gradient_for_glob_conv_newton_solve |
DoubleVector storing the gradient for the globally convergent Newton method. | |
Dense LU decomposition-based solve of linear system assembled via finite differencing of the residuals Vector. Even more inefficient than DenseLU but excellent sanity check!
Definition at line 433 of file linear_solver.h.
|
inline |
Constructor: empty.
Definition at line 437 of file linear_solver.h.
|
inlinevirtual |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system.
Reimplemented from oomph::DenseLU.
Definition at line 452 of file linear_solver.h.
References oomph::DenseLU::solve().
|
inlinevirtual |
Linear-algebra-type solver: Takes pointer to a matrix and rhs vector and returns the solution of the linear system Call the broken base-class version. If you want this, please implement it.
Reimplemented from oomph::DenseLU.
Definition at line 463 of file linear_solver.h.
References oomph::LinearSolver::solve().
|
virtual |
Solver: Takes pointer to problem and returns the results Vector which contains the solution of the linear system defined by the problem's residual Vector (Jacobian computed by FD approx.)
Solver: Takes pointer to problem and returns the results Vector which contains the solution of the linear system defined by the problem's residual Vector. (Jacobian assembled by FD).
Reimplemented from oomph::DenseLU.
Definition at line 579 of file linear_solver.cc.
References oomph::TimingHelpers::convert_secs_to_formatted_string(), oomph::LinearSolver::Doc_time, oomph::Problem::get_fd_jacobian(), oomph::DenseLU::Jacobian_setup_time, oomph::Problem::ndof(), oomph::oomph_info, oomph::DenseLU::Sign_of_determinant_of_matrix, oomph::Problem::sign_of_jacobian(), and solve().
Referenced by solve().