Public Types | Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
oomph::MGSolver< DIM > Class Template Reference

/////////////////////////////////////////////////////// /////////////////////////////////////////////////////// More...

#include <geometric_multigrid.h>

+ Inheritance diagram for oomph::MGSolver< DIM >:

Public Types

typedef Smoother *(* PreSmootherFactoryFctPt) ()
 typedef for a function that returns a pointer to an object of the class Smoother to be used as the pre-smoother More...
 
typedef Smoother *(* PostSmootherFactoryFctPt) ()
 typedef for a function that returns a pointer to an object of the class Smoother to be used as the post-smoother More...
 

Public Member Functions

void set_pre_smoother_factory_function (PreSmootherFactoryFctPt pre_smoother_fn)
 Access function to set the pre-smoother creation function. More...
 
void set_post_smoother_factory_function (PostSmootherFactoryFctPt post_smoother_fn)
 Access function to set the post-smoother creation function. More...
 
 MGSolver (MGProblem *mg_problem_pt)
 Constructor: Set up default values for number of V-cycles and pre- and post-smoothing steps. More...
 
 ~MGSolver ()
 Delete any dynamically allocated data. More...
 
void clean_up_memory ()
 Clean up anything that needs to be cleaned up. More...
 
void set_self_test_vector ()
 Makes a vector which will be used in the self-test. Is currently set to make the entries of the vector represent a plane wave propagating at an angle of 45 degrees. More...
 
void self_test ()
 Makes a vector, restricts it down the levels of the hierarchy and documents it at each level. After this is done the vector is interpolated up the levels of the hierarchy with the output being documented at each level. More...
 
void restriction_self_test ()
 Make a self-test to make sure that the interpolation matrices are doing the same thing to restrict the vectors down through the heirachy. More...
 
void interpolation_self_test ()
 Make a self-test to make sure that the interpolation matrices are doing the same thing to interpolate the vectors up. More...
 
void plot (const unsigned &hierarchy_level, const DoubleVector &input_vector, const std::string &filename)
 Given a level in the hierarchy, an input vector and a filename this function will document the given vector according to the structure of the mesh on the given level. More...
 
void disable_v_cycle_output ()
 Disable all output from mg_solve apart from the number of V-cycles used to solve the problem. More...
 
void disable_output ()
 Suppress anything that can be suppressed, i.e. any timings. Things like mesh adaptation can not however be silenced using this. More...
 
void enable_v_cycle_output ()
 Enable the output of the V-cycle timings and other output. More...
 
void enable_doc_everything ()
 Enable the output from anything that could have been suppressed. More...
 
void enable_output ()
 Enable the output from anything that could have been suppressed. More...
 
void disable_smoother_and_superlu_doc_time ()
 Suppress the output of both smoothers and SuperLU. More...
 
unsigned & npost_smooth ()
 Return the number of post-smoothing iterations (lvalue) More...
 
unsigned & npre_smooth ()
 Return the number of pre-smoothing iterations (lvalue) More...
 
void pre_smooth (const unsigned &level)
 Pre-smoother: Perform 'max_iter' smoothing steps on the linear system Ax=b with current RHS vector, b, starting with current solution vector, x. Return the residual vector r=b-Ax. Uses the default smoother (set in the MGProblem constructor) which can be overloaded for a specific problem. More...
 
void post_smooth (const unsigned &level)
 Post-smoother: Perform max_iter smoothing steps on the linear system Ax=b with current RHS vector, b, starting with current solution vector, x. Uses the default smoother (set in the MGProblem constructor) which can be overloaded for specific problem. More...
 
double residual_norm (const unsigned &level)
 Return norm of residual r=b-Ax and the residual vector itself on the level-th level. More...
 
void direct_solve ()
 Call the direct solver (SuperLU) to solve the problem exactly. More...
 
void interpolation_matrix_set (const unsigned &level, double *value, int *col_index, int *row_st, unsigned &ncol, unsigned &nnz)
 Builds a CRDoubleMatrix that is used to interpolate the residual between levels. The transpose can be used as the full weighting restriction. More...
 
void interpolation_matrix_set (const unsigned &level, Vector< double > &value, Vector< int > &col_index, Vector< int > &row_st, unsigned &ncol, unsigned &nrow)
 Builds a CRDoubleMatrix that is used to interpolate the residual between levels. The transpose can be used as the full weighting restriction. More...
 
void set_restriction_matrices_as_interpolation_transposes ()
 Builds a CRDoubleMatrix on each level that is used to restrict the residual between levels. The transpose can be used as the interpolation matrix. More...
 
void restrict_residual (const unsigned &level)
 Restrict residual (computed on level-th MG level) to the next coarser mesh and stick it into the coarse mesh RHS vector. More...
 
void interpolate_and_correct (const unsigned &level)
 Interpolate solution at current level onto next finer mesh and correct the solution x at that level. More...
 
void level_up_local_coord_of_node (const int &son_type, Vector< double > &s)
 Given the son_type of an element and a local node number j in that element with nnode_1d nodes per coordinate direction, return the local coordinate s in its father element. Needed in the setup of the interpolation matrices. More...
 
void setup_interpolation_matrices ()
 Setup the interpolation matrix on each level. More...
 
void setup_interpolation_matrices_unstructured ()
 Setup the interpolation matrix on each level (used for unstructured meshes) More...
 
void setup_transfer_matrices ()
 Setup the transfer matrices on each level. More...
 
void full_setup ()
 Do a full setup (assumes everything will be setup around the MGProblem pointer given in the constructor) More...
 
void solve (Problem *const &problem_pt, DoubleVector &result)
 Virtual function in the base class that needs to be implemented later but for now just leave it empty. More...
 
unsigned iterations () const
 Number of iterations. More...
 
unsigned & max_iter ()
 Number of iterations. More...
 
void level_up_local_coord_of_node (const int &son_type, Vector< double > &s)
 Given the son type of the element and the local coordinate s of a given node in the son element, return the local coordinate s in its father element. 2D case. More...
 
void level_up_local_coord_of_node (const int &son_type, Vector< double > &s)
 Given the son type of the element and the local coordinate s of a given node in the son element, return the local coordinate s in its father element. 3D case. More...
 
- Public Member Functions inherited from oomph::IterativeLinearSolver
 IterativeLinearSolver ()
 Constructor: Set (default) trivial preconditioner and set defaults for tolerance and max. number of iterations. More...
 
 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="")
 Write convergence history into file with specified filename (automatically switches on doc). Optional second argument is a string that can be used (as a zone title) to identify what case we're running (e.g. what combination of linear solver and preconditioner or parameter values are used). More...
 
void close_convergence_history_file_stream ()
 Close convergence history output stream. More...
 
double jacobian_setup_time () const
 returns the time taken to assemble the jacobian matrix and residual vector More...
 
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 ()
 Enables the iterative solver be used as preconditioner (when calling the solve method it bypass the setup solver method – currently only used by Trilinos solver —) More...
 
void disable_iterative_solver_as_preconditioner ()
 Disables the iterative solver be used as preconditioner (when calling the solve method it bypass the setup solver method – currently only used by Trilinos solver —) More...
 
- Public Member Functions inherited from oomph::LinearSolver
 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 ()
 Enable resolve (i.e. store matrix and/or LU decomposition, say) Virtual so it can be overloaded to perform additional tasks. More...
 
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. More...
 
virtual 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. More...
 
virtual 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. More...
 
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). More...
 
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. More...
 
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. More...
 
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) More...
 
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) More...
 
virtual void enable_computation_of_gradient ()
 function to enable the computation of the gradient required for the globally convergent Newton method More...
 
void disable_computation_of_gradient ()
 function to disable the computation of the gradient required for the globally convergent Newton method More...
 
void reset_gradient ()
 function to reset the size of the gradient before each Newton solve More...
 
void get_gradient (DoubleVector &gradient)
 function to access the gradient, provided it has been computed More...
 
- Public Member Functions inherited from oomph::DistributableLinearAlgebraObject
 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...
 
LinearAlgebraDistributiondistribution_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
 if the communicator_pt is null then the distribution is not setup then false is returned, otherwise return true More...
 
void build_distribution (const LinearAlgebraDistribution *const dist_pt)
 setup the distribution of this distributable linear algebra object More...
 
void build_distribution (const LinearAlgebraDistribution &dist)
 setup the distribution of this distributable linear algebra object More...
 

Protected Member Functions

void mg_solve (DoubleVector &result)
 Do the actual solve – this is called through the pure virtual solve function in the LinearSolver base class. The function is stored as protected to allow the MGPreconditioner derived class to use the solver. More...
 
void modify_restriction_matrices ()
 Normalise the rows of the restriction matrices to avoid amplifications when projecting to the coarser level. More...
 
- Protected Member Functions inherited from oomph::DistributableLinearAlgebraObject
void clear_distribution ()
 clear the distribution of this distributable linear algebra object More...
 

Protected Attributes

unsigned Nvcycle
 Maximum number of V-cycles (this is set as a protected variable so. More...
 
MGProblemMg_problem_pt
 Pointer to the MG problem (deep copy). This is protected to provide access to the MG preconditioner. More...
 
Vector< DoubleVectorRhs_mg_vectors_storage
 Vector to store the RHS vectors (Rhs_mg). This is protected to allow the multigrid preconditioner to assign the RHS vector during preconditioner_solve() More...
 
bool Suppress_v_cycle_output
 Indicates whether or not the V-cycle output should be suppressed. Needs to be protected member data for the multigrid preconditioner to know whether or not to output information with each preconditioning step. More...
 
bool Suppress_all_output
 If this is set to true then all output from the solver is suppressed. This is protected member data so that the multigrid preconditioner knows whether or not to restore the stream pointer. More...
 
std::ostream * Stream_pt
 Pointer to the output stream – defaults to std::cout. This is protected member data to allow the preconditioner to restore normal output if everything was chosen to be suppressed by the user. More...
 
- Protected Attributes inherited from oomph::IterativeLinearSolver
bool Doc_convergence_history
 Flag indicating if the convergence history is to be documented. More...
 
std::ofstream Output_file_stream
 Output file stream for convergence history. More...
 
double Tolerance
 Convergence tolerance. More...
 
unsigned Max_iter
 Maximum number of iterations. More...
 
PreconditionerPreconditioner_pt
 Pointer to the preconditioner. More...
 
double Jacobian_setup_time
 Jacobian setup time. More...
 
double Solution_time
 linear solver solution time More...
 
double Preconditioner_setup_time
 Preconditioner setup time. More...
 
bool Setup_preconditioner_before_solve
 indicates whether the preconditioner should be setup before solve. Default = true; More...
 
bool Throw_error_after_max_iter
 Should we throw an error instead of just returning when we hit the max iterations? More...
 
bool Use_iterative_solver_as_preconditioner
 Use the iterative solver as preconditioner. More...
 
bool First_time_solve_when_used_as_preconditioner
 When the iterative solver is used a preconditioner then we call the setup of solver method only once (the first time the solve method is called) More...
 
- Protected Attributes inherited from oomph::LinearSolver
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. More...
 
bool Doc_time
 Boolean flag that indicates whether the time taken. More...
 
bool Compute_gradient
 flag that indicates whether the gradient required for the globally convergent Newton method should be computed or not More...
 
bool Gradient_has_been_computed
 flag that indicates whether the gradient was computed or not More...
 
DoubleVector Gradient_for_glob_conv_newton_solve
 DoubleVector storing the gradient for the globally convergent Newton method. More...
 

Private Member Functions

void setup_mg_hierarchy ()
 Function to set up the hierachy of levels. Creates a vector of pointers to each MG level. More...
 
void setup_mg_structures ()
 Function to set up the hierachy of levels. Creates a vector of pointers to each MG level. More...
 
void setup_smoothers ()
 Function to set up all of the smoothers once the system matrices have been set up. More...
 

Private Attributes

PreSmootherFactoryFctPt Pre_smoother_factory_function_pt
 Function to create pre-smoothers. More...
 
PostSmootherFactoryFctPt Post_smoother_factory_function_pt
 Function to create post-smoothers. More...
 
unsigned Nlevel
 The number of levels in the multigrid heirachy. More...
 
Vector< MGProblem * > Mg_hierarchy
 Vector containing pointers to problems in hierarchy. More...
 
Vector< CRDoubleMatrix * > Mg_matrices_storage_pt
 Vector to store the system matrices. More...
 
Vector< CRDoubleMatrix * > Interpolation_matrices_storage_pt
 Vector to store the interpolation matrices. More...
 
Vector< CRDoubleMatrix * > Restriction_matrices_storage_pt
 Vector to store the restriction matrices. More...
 
Vector< DoubleVectorX_mg_vectors_storage
 Vector to store the solution vectors (X_mg) More...
 
Vector< DoubleVectorResidual_mg_vectors_storage
 Vector to store the residual vectors. More...
 
Vector< DoubleVectorInterpolation_self_test_vectors_storage
 Vector to store the result of interpolation on each level (only required if the user wishes to document the output of interpolation and restriction on each level) More...
 
Vector< DoubleVectorRestriction_self_test_vectors_storage
 Vector to store the result of restriction on each level (only required if the user wishes to document the output of interpolation and restriction on each level) More...
 
Vector< Smoother * > Pre_smoothers_storage_pt
 Vector to store the pre-smoothers. More...
 
Vector< Smoother * > Post_smoothers_storage_pt
 Vector to store the post-smoothers. More...
 
unsigned Npre_smooth
 Number of pre-smoothing steps. More...
 
unsigned Npost_smooth
 Number of post-smoothing steps. More...
 
bool Doc_everything
 If this is set to true we document everything. In addition to outputting the information of the setup timings and V-cycle data we document the refinement and unrefinement patterns given by the transfer operators which is done by keeping the coarser MG problem pointers alive. More...
 
bool Has_been_setup
 Boolean variable to indicate whether or not the solver has been setup. More...
 
bool Has_been_solved
 Boolean variable to indicate whether or not the problem was successfully solved. More...
 
unsigned V_cycle_counter
 Pointer to counter for V-cycles. More...
 

Additional Inherited Members

- Static Protected Attributes inherited from oomph::IterativeLinearSolver
static IdentityPreconditioner Default_preconditioner
 Default preconditioner: The base class for preconditioners is a fully functional (if trivial!) preconditioner. More...
 

Detailed Description

template<unsigned DIM>
class oomph::MGSolver< DIM >

/////////////////////////////////////////////////////// ///////////////////////////////////////////////////////

Definition at line 89 of file geometric_multigrid.h.

Member Typedef Documentation

◆ PostSmootherFactoryFctPt

template<unsigned DIM>
typedef Smoother*(* oomph::MGSolver< DIM >::PostSmootherFactoryFctPt) ()

typedef for a function that returns a pointer to an object of the class Smoother to be used as the post-smoother

Definition at line 98 of file geometric_multigrid.h.

◆ PreSmootherFactoryFctPt

template<unsigned DIM>
typedef Smoother*(* oomph::MGSolver< DIM >::PreSmootherFactoryFctPt) ()

typedef for a function that returns a pointer to an object of the class Smoother to be used as the pre-smoother

Definition at line 94 of file geometric_multigrid.h.

Constructor & Destructor Documentation

◆ MGSolver()

template<unsigned DIM>
oomph::MGSolver< DIM >::MGSolver ( MGProblem mg_problem_pt)
inline

Constructor: Set up default values for number of V-cycles and pre- and post-smoothing steps.

Definition at line 118 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Mg_hierarchy, oomph::MGSolver< DIM >::Mg_problem_pt, and oomph::IterativeLinearSolver::Tolerance.

◆ ~MGSolver()

template<unsigned DIM>
oomph::MGSolver< DIM >::~MGSolver ( )
inline

Delete any dynamically allocated data.

Definition at line 141 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::clean_up_memory().

Member Function Documentation

◆ clean_up_memory()

template<unsigned DIM>
void oomph::MGSolver< DIM >::clean_up_memory ( )
inlinevirtual

◆ direct_solve()

template<unsigned DIM>
void oomph::MGSolver< DIM >::direct_solve ( )
inline

◆ disable_output()

template<unsigned DIM>
void oomph::MGSolver< DIM >::disable_output ( )
inline

Suppress anything that can be suppressed, i.e. any timings. Things like mesh adaptation can not however be silenced using this.

Definition at line 256 of file geometric_multigrid.h.

References oomph::LinearSolver::Doc_time, oomph::oomph_info, oomph::oomph_nullstream, oomph::MGSolver< DIM >::Stream_pt, oomph::OomphInfo::stream_pt(), oomph::MGSolver< DIM >::Suppress_all_output, and oomph::MGSolver< DIM >::Suppress_v_cycle_output.

◆ disable_smoother_and_superlu_doc_time()

template<unsigned DIM>
void oomph::MGSolver< DIM >::disable_smoother_and_superlu_doc_time ( )
inline

◆ disable_v_cycle_output()

template<unsigned DIM>
void oomph::MGSolver< DIM >::disable_v_cycle_output ( )
inline

Disable all output from mg_solve apart from the number of V-cycles used to solve the problem.

Definition at line 245 of file geometric_multigrid.h.

References oomph::LinearSolver::Doc_time, and oomph::MGSolver< DIM >::Suppress_v_cycle_output.

◆ enable_doc_everything()

template<unsigned DIM>
void oomph::MGSolver< DIM >::enable_doc_everything ( )
inline

Enable the output from anything that could have been suppressed.

Definition at line 286 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Doc_everything.

◆ enable_output()

template<unsigned DIM>
void oomph::MGSolver< DIM >::enable_output ( )
inline

Enable the output from anything that could have been suppressed.

Definition at line 295 of file geometric_multigrid.h.

References oomph::LinearSolver::Doc_time, oomph::MGSolver< DIM >::Suppress_all_output, and oomph::MGSolver< DIM >::Suppress_v_cycle_output.

◆ enable_v_cycle_output()

template<unsigned DIM>
void oomph::MGSolver< DIM >::enable_v_cycle_output ( )
inline

Enable the output of the V-cycle timings and other output.

Definition at line 276 of file geometric_multigrid.h.

References oomph::LinearSolver::Doc_time, and oomph::MGSolver< DIM >::Suppress_v_cycle_output.

◆ full_setup()

template<unsigned DIM>
void oomph::MGSolver< DIM >::full_setup

Do a full setup (assumes everything will be setup around the MGProblem pointer given in the constructor)

Runs a full setup of the MG solver.

Definition at line 826 of file geometric_multigrid.h.

References oomph::TerminateHelper::clean_up_memory(), oomph::FiniteElement::dim(), i, oomph::oomph_info, oomph::Global_string_for_annotation::string(), and oomph::TimingHelpers::timer().

Referenced by oomph::MGPreconditioner< DIM >::setup(), and oomph::MGSolver< DIM >::solve().

◆ interpolate_and_correct()

template<unsigned DIM>
void oomph::MGSolver< DIM >::interpolate_and_correct ( const unsigned &  level)

Interpolate solution at current level onto next finer mesh and correct the solution x at that level.

Definition at line 2043 of file geometric_multigrid.h.

◆ interpolation_matrix_set() [1/2]

template<unsigned DIM>
void oomph::MGSolver< DIM >::interpolation_matrix_set ( const unsigned &  level,
double *  value,
int *  col_index,
int *  row_st,
unsigned &  ncol,
unsigned &  nnz 
)
inline

Builds a CRDoubleMatrix that is used to interpolate the residual between levels. The transpose can be used as the full weighting restriction.

Definition at line 402 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Interpolation_matrices_storage_pt.

◆ interpolation_matrix_set() [2/2]

template<unsigned DIM>
void oomph::MGSolver< DIM >::interpolation_matrix_set ( const unsigned &  level,
Vector< double > &  value,
Vector< int > &  col_index,
Vector< int > &  row_st,
unsigned &  ncol,
unsigned &  nrow 
)
inline

Builds a CRDoubleMatrix that is used to interpolate the residual between levels. The transpose can be used as the full weighting restriction.

Definition at line 420 of file geometric_multigrid.h.

References oomph::LinearAlgebraDistribution::communicator_pt(), oomph::MGSolver< DIM >::Interpolation_matrices_storage_pt, oomph::MGSolver< DIM >::Mg_hierarchy, oomph::DistributableLinearAlgebraObject::nrow(), and oomph::Global_string_for_annotation::string().

◆ interpolation_self_test()

template<unsigned DIM>
void oomph::MGSolver< DIM >::interpolation_self_test

Make a self-test to make sure that the interpolation matrices are doing the same thing to interpolate the vectors up.

Function which implements a self-test to interpolate a vector up all of levels of the MG hierarchy and outputs the restricted vectors to file.

Definition at line 2344 of file geometric_multigrid.h.

References oomph::Global_string_for_annotation::string().

◆ iterations()

template<unsigned DIM>
unsigned oomph::MGSolver< DIM >::iterations ( ) const
inlinevirtual

Number of iterations.

Implements oomph::IterativeLinearSolver.

Definition at line 595 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::V_cycle_counter.

◆ level_up_local_coord_of_node() [1/3]

void oomph::MGSolver< 2 >::level_up_local_coord_of_node ( const int &  son_type,
Vector< double > &  s 
)

Given the son type of the element and the local coordinate s of a given node in the son element, return the local coordinate s in its father element. 2D case.

Definition at line 37 of file geometric_multigrid.cc.

References oomph::QuadTreeNames::NE, oomph::QuadTreeNames::NW, oomph::Tree::OMEGA, s, oomph::QuadTreeNames::SE, and oomph::QuadTreeNames::SW.

◆ level_up_local_coord_of_node() [2/3]

void oomph::MGSolver< 3 >::level_up_local_coord_of_node ( const int &  son_type,
Vector< double > &  s 
)

Given the son type of the element and the local coordinate s of a given node in the son element, return the local coordinate s in its father element. 3D case.

Definition at line 91 of file geometric_multigrid.cc.

References oomph::OcTreeNames::LDB, oomph::OcTreeNames::LDF, oomph::OcTreeNames::LUB, oomph::OcTreeNames::LUF, oomph::Tree::OMEGA, oomph::OcTreeNames::RDB, oomph::OcTreeNames::RDF, oomph::OcTreeNames::RUB, oomph::OcTreeNames::RUF, and s.

◆ level_up_local_coord_of_node() [3/3]

template<unsigned DIM>
void oomph::MGSolver< DIM >::level_up_local_coord_of_node ( const int &  son_type,
Vector< double > &  s 
)

Given the son_type of an element and a local node number j in that element with nnode_1d nodes per coordinate direction, return the local coordinate s in its father element. Needed in the setup of the interpolation matrices.

◆ max_iter()

template<unsigned DIM>
unsigned& oomph::MGSolver< DIM >::max_iter ( )
inline

Number of iterations.

Definition at line 602 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Nvcycle.

◆ mg_solve()

template<unsigned DIM>
void oomph::MGSolver< DIM >::mg_solve ( DoubleVector result)
protected

Do the actual solve – this is called through the pure virtual solve function in the LinearSolver base class. The function is stored as protected to allow the MGPreconditioner derived class to use the solver.

Linear solver.

Definition at line 2525 of file geometric_multigrid.h.

References i, oomph::oomph_info, and oomph::TimingHelpers::timer().

Referenced by oomph::MGPreconditioner< DIM >::preconditioner_solve(), and oomph::MGSolver< DIM >::solve().

◆ modify_restriction_matrices()

template<unsigned DIM>
void oomph::MGSolver< DIM >::modify_restriction_matrices
protected

Normalise the rows of the restriction matrices to avoid amplifications when projecting to the coarser level.

Modify the restriction matrices.

Definition at line 2070 of file geometric_multigrid.h.

References i, oomph::CRDoubleMatrix::nrow(), oomph::CRDoubleMatrix::row_start(), and oomph::CRDoubleMatrix::value().

◆ npost_smooth()

template<unsigned DIM>
unsigned& oomph::MGSolver< DIM >::npost_smooth ( )
inline

Return the number of post-smoothing iterations (lvalue)

Definition at line 328 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Npost_smooth.

◆ npre_smooth()

template<unsigned DIM>
unsigned& oomph::MGSolver< DIM >::npre_smooth ( )
inline

Return the number of pre-smoothing iterations (lvalue)

Definition at line 336 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Npre_smooth.

◆ plot()

template<unsigned DIM>
void oomph::MGSolver< DIM >::plot ( const unsigned &  hierarchy_level,
const DoubleVector input_vector,
const std::string &  filename 
)

Given a level in the hierarchy, an input vector and a filename this function will document the given vector according to the structure of the mesh on the given level.

Plots the input vector (assuming we're dealing with scalar nodal data, otherwise I don't know how to implement this...)

Definition at line 2379 of file geometric_multigrid.h.

References e, oomph::Data::eqn_number(), oomph::Mesh::finite_element_pt(), oomph::Node::hanging_pt(), i, oomph::Node::is_hanging(), oomph::Node::is_on_boundary(), oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::HangInfo::nmaster(), oomph::FiniteElement::nnode(), oomph::FiniteElement::nnode_1d(), oomph::FiniteElement::node_pt(), oomph::Mesh::node_pt(), oomph::Data::nvalue(), oomph::FiniteElement::tecplot_zone_string(), and oomph::Node::x().

◆ post_smooth()

template<unsigned DIM>
void oomph::MGSolver< DIM >::post_smooth ( const unsigned &  level)
inline

Post-smoother: Perform max_iter smoothing steps on the linear system Ax=b with current RHS vector, b, starting with current solution vector, x. Uses the default smoother (set in the MGProblem constructor) which can be overloaded for specific problem.

Definition at line 366 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Post_smoothers_storage_pt, oomph::MGSolver< DIM >::Rhs_mg_vectors_storage, and oomph::MGSolver< DIM >::X_mg_vectors_storage.

◆ pre_smooth()

template<unsigned DIM>
void oomph::MGSolver< DIM >::pre_smooth ( const unsigned &  level)
inline

Pre-smoother: Perform 'max_iter' smoothing steps on the linear system Ax=b with current RHS vector, b, starting with current solution vector, x. Return the residual vector r=b-Ax. Uses the default smoother (set in the MGProblem constructor) which can be overloaded for a specific problem.

Definition at line 348 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Mg_matrices_storage_pt, oomph::MGSolver< DIM >::Pre_smoothers_storage_pt, oomph::MGSolver< DIM >::Residual_mg_vectors_storage, oomph::MGSolver< DIM >::Rhs_mg_vectors_storage, and oomph::MGSolver< DIM >::X_mg_vectors_storage.

◆ residual_norm()

template<unsigned DIM>
double oomph::MGSolver< DIM >::residual_norm ( const unsigned &  level)
inline

◆ restrict_residual()

template<unsigned DIM>
void oomph::MGSolver< DIM >::restrict_residual ( const unsigned &  level)

Restrict residual (computed on level-th MG level) to the next coarser mesh and stick it into the coarse mesh RHS vector.

Restrict residual (computed on current MG level) to next coarser mesh and stick it into the coarse mesh RHS vector using the restriction matrix (if restrict_flag=1) or the transpose of the interpolation matrix (if restrict_flag=2)

Definition at line 2020 of file geometric_multigrid.h.

◆ restriction_self_test()

template<unsigned DIM>
void oomph::MGSolver< DIM >::restriction_self_test

Make a self-test to make sure that the interpolation matrices are doing the same thing to restrict the vectors down through the heirachy.

Function which implements a self-test to restrict the residual vector down all of the coarser grids and output the restricted vectors to file.

Definition at line 2307 of file geometric_multigrid.h.

References oomph::Global_string_for_annotation::string().

◆ self_test()

template<unsigned DIM>
void oomph::MGSolver< DIM >::self_test

Makes a vector, restricts it down the levels of the hierarchy and documents it at each level. After this is done the vector is interpolated up the levels of the hierarchy with the output being documented at each level.

Definition at line 2135 of file geometric_multigrid.h.

References oomph::LinearAlgebraDistribution::communicator_pt(), oomph::oomph_info, oomph::Global_string_for_annotation::string(), and oomph::TimingHelpers::timer().

◆ set_post_smoother_factory_function()

template<unsigned DIM>
void oomph::MGSolver< DIM >::set_post_smoother_factory_function ( PostSmootherFactoryFctPt  post_smoother_fn)
inline

Access function to set the post-smoother creation function.

Definition at line 109 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Post_smoother_factory_function_pt.

◆ set_pre_smoother_factory_function()

template<unsigned DIM>
void oomph::MGSolver< DIM >::set_pre_smoother_factory_function ( PreSmootherFactoryFctPt  pre_smoother_fn)
inline

Access function to set the pre-smoother creation function.

Definition at line 101 of file geometric_multigrid.h.

References oomph::MGSolver< DIM >::Pre_smoother_factory_function_pt.

◆ set_restriction_matrices_as_interpolation_transposes()

template<unsigned DIM>
void oomph::MGSolver< DIM >::set_restriction_matrices_as_interpolation_transposes ( )
inline

Builds a CRDoubleMatrix on each level that is used to restrict the residual between levels. The transpose can be used as the interpolation matrix.

Definition at line 465 of file geometric_multigrid.h.

References i, oomph::MGSolver< DIM >::Interpolation_matrices_storage_pt, oomph::MGSolver< DIM >::Nlevel, and oomph::MGSolver< DIM >::Restriction_matrices_storage_pt.

◆ set_self_test_vector()

template<unsigned DIM>
void oomph::MGSolver< DIM >::set_self_test_vector

Makes a vector which will be used in the self-test. Is currently set to make the entries of the vector represent a plane wave propagating at an angle of 45 degrees.

Sets the initial vector to be used in the restriction and interpolation self-tests.

Definition at line 2219 of file geometric_multigrid.h.

References e, oomph::Data::eqn_number(), oomph::Mesh::finite_element_pt(), i, oomph::Node::ndim(), oomph::Mesh::nelement(), oomph::FiniteElement::nnode(), oomph::FiniteElement::node_pt(), oomph::Mesh::node_pt(), oomph::Data::nvalue(), oomph::MathematicalConstants::Pi, and oomph::Node::x().

◆ setup_interpolation_matrices()

template<unsigned DIM>
void oomph::MGSolver< DIM >::setup_interpolation_matrices

Setup the interpolation matrix on each level.

Setup the interpolation matrices.

Definition at line 1510 of file geometric_multigrid.h.

References oomph::Data::eqn_number(), oomph::Mesh::finite_element_pt(), i, oomph::HangInfo::master_node_pt(), oomph::HangInfo::master_weight(), oomph::Mesh::nelement(), oomph::HangInfo::nmaster(), oomph::Tree::OMEGA, and s.

◆ setup_interpolation_matrices_unstructured()

template<unsigned DIM>
void oomph::MGSolver< DIM >::setup_interpolation_matrices_unstructured

◆ setup_mg_hierarchy()

template<unsigned DIM>
void oomph::MGSolver< DIM >::setup_mg_hierarchy
private

◆ setup_mg_structures()

template<unsigned DIM>
void oomph::MGSolver< DIM >::setup_mg_structures
private

Function to set up the hierachy of levels. Creates a vector of pointers to each MG level.

Set up the MG hierarchy structures.

Definition at line 1180 of file geometric_multigrid.h.

References oomph::LinearAlgebraDistribution::clear(), oomph::LinearAlgebraDistribution::communicator_pt(), i, oomph::oomph_info, oomph::Global_string_for_annotation::string(), and oomph::TimingHelpers::timer().

◆ setup_smoothers()

template<unsigned DIM>
void oomph::MGSolver< DIM >::setup_smoothers
private

Function to set up all of the smoothers once the system matrices have been set up.

Set up the smoothers on all levels.

Definition at line 1360 of file geometric_multigrid.h.

References oomph::LinearAlgebraDistribution::communicator_pt(), i, oomph::oomph_info, oomph::Global_string_for_annotation::string(), oomph::TimingHelpers::timer(), and oomph::IterativeLinearSolver::tolerance().

◆ setup_transfer_matrices()

template<unsigned DIM>
void oomph::MGSolver< DIM >::setup_transfer_matrices

Setup the transfer matrices on each level.

Set up the transfer matrices. Both the pure injection and full weighting method have been implemented here but it is highly recommended that full weighting is used in general. In both methods the transpose of the transfer matrix is used to transfer a vector back.

Definition at line 1117 of file geometric_multigrid.h.

References oomph::oomph_info, and oomph::TimingHelpers::timer().

◆ solve()

template<unsigned DIM>
void oomph::MGSolver< DIM >::solve ( Problem *const &  problem_pt,
DoubleVector result 
)
inlinevirtual

Member Data Documentation

◆ Doc_everything

template<unsigned DIM>
bool oomph::MGSolver< DIM >::Doc_everything
private

If this is set to true we document everything. In addition to outputting the information of the setup timings and V-cycle data we document the refinement and unrefinement patterns given by the transfer operators which is done by keeping the coarser MG problem pointers alive.

Definition at line 717 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::clean_up_memory(), and oomph::MGSolver< DIM >::enable_doc_everything().

◆ Has_been_setup

template<unsigned DIM>
bool oomph::MGSolver< DIM >::Has_been_setup
private

Boolean variable to indicate whether or not the solver has been setup.

Definition at line 720 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::clean_up_memory().

◆ Has_been_solved

template<unsigned DIM>
bool oomph::MGSolver< DIM >::Has_been_solved
private

Boolean variable to indicate whether or not the problem was successfully solved.

Definition at line 724 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::solve().

◆ Interpolation_matrices_storage_pt

template<unsigned DIM>
Vector<CRDoubleMatrix*> oomph::MGSolver< DIM >::Interpolation_matrices_storage_pt
private

◆ Interpolation_self_test_vectors_storage

template<unsigned DIM>
Vector<DoubleVector> oomph::MGSolver< DIM >::Interpolation_self_test_vectors_storage
private

Vector to store the result of interpolation on each level (only required if the user wishes to document the output of interpolation and restriction on each level)

Definition at line 693 of file geometric_multigrid.h.

◆ Mg_hierarchy

template<unsigned DIM>
Vector<MGProblem*> oomph::MGSolver< DIM >::Mg_hierarchy
private

◆ Mg_matrices_storage_pt

template<unsigned DIM>
Vector<CRDoubleMatrix*> oomph::MGSolver< DIM >::Mg_matrices_storage_pt
private

◆ Mg_problem_pt

template<unsigned DIM>
MGProblem* oomph::MGSolver< DIM >::Mg_problem_pt
protected

Pointer to the MG problem (deep copy). This is protected to provide access to the MG preconditioner.

Definition at line 626 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::MGSolver(), oomph::MGPreconditioner< DIM >::preconditioner_solve(), and oomph::MGSolver< DIM >::solve().

◆ Nlevel

template<unsigned DIM>
unsigned oomph::MGSolver< DIM >::Nlevel
private

◆ Npost_smooth

template<unsigned DIM>
unsigned oomph::MGSolver< DIM >::Npost_smooth
private

Number of post-smoothing steps.

Definition at line 710 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::npost_smooth().

◆ Npre_smooth

template<unsigned DIM>
unsigned oomph::MGSolver< DIM >::Npre_smooth
private

Number of pre-smoothing steps.

Definition at line 707 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::npre_smooth().

◆ Nvcycle

template<unsigned DIM>
unsigned oomph::MGSolver< DIM >::Nvcycle
protected

Maximum number of V-cycles (this is set as a protected variable so.

Definition at line 622 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::max_iter(), and oomph::MGPreconditioner< DIM >::MGPreconditioner().

◆ Post_smoother_factory_function_pt

template<unsigned DIM>
PostSmootherFactoryFctPt oomph::MGSolver< DIM >::Post_smoother_factory_function_pt
private

Function to create post-smoothers.

Definition at line 655 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::set_post_smoother_factory_function().

◆ Post_smoothers_storage_pt

template<unsigned DIM>
Vector<Smoother*> oomph::MGSolver< DIM >::Post_smoothers_storage_pt
private

◆ Pre_smoother_factory_function_pt

template<unsigned DIM>
PreSmootherFactoryFctPt oomph::MGSolver< DIM >::Pre_smoother_factory_function_pt
private

Function to create pre-smoothers.

Definition at line 652 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::set_pre_smoother_factory_function().

◆ Pre_smoothers_storage_pt

template<unsigned DIM>
Vector<Smoother*> oomph::MGSolver< DIM >::Pre_smoothers_storage_pt
private

◆ Residual_mg_vectors_storage

template<unsigned DIM>
Vector<DoubleVector> oomph::MGSolver< DIM >::Residual_mg_vectors_storage
private

Vector to store the residual vectors.

Definition at line 688 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::pre_smooth(), and oomph::MGSolver< DIM >::residual_norm().

◆ Restriction_matrices_storage_pt

template<unsigned DIM>
Vector<CRDoubleMatrix*> oomph::MGSolver< DIM >::Restriction_matrices_storage_pt
private

◆ Restriction_self_test_vectors_storage

template<unsigned DIM>
Vector<DoubleVector> oomph::MGSolver< DIM >::Restriction_self_test_vectors_storage
private

Vector to store the result of restriction on each level (only required if the user wishes to document the output of interpolation and restriction on each level)

Definition at line 698 of file geometric_multigrid.h.

◆ Rhs_mg_vectors_storage

template<unsigned DIM>
Vector<DoubleVector> oomph::MGSolver< DIM >::Rhs_mg_vectors_storage
protected

Vector to store the RHS vectors (Rhs_mg). This is protected to allow the multigrid preconditioner to assign the RHS vector during preconditioner_solve()

Definition at line 631 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::direct_solve(), oomph::MGSolver< DIM >::post_smooth(), oomph::MGSolver< DIM >::pre_smooth(), oomph::MGPreconditioner< DIM >::preconditioner_solve(), and oomph::MGSolver< DIM >::residual_norm().

◆ Stream_pt

template<unsigned DIM>
std::ostream* oomph::MGSolver< DIM >::Stream_pt
protected

Pointer to the output stream – defaults to std::cout. This is protected member data to allow the preconditioner to restore normal output if everything was chosen to be suppressed by the user.

Definition at line 648 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::disable_output(), oomph::MGPreconditioner< DIM >::preconditioner_solve(), oomph::MGPreconditioner< DIM >::setup(), and oomph::MGSolver< DIM >::solve().

◆ Suppress_all_output

template<unsigned DIM>
bool oomph::MGSolver< DIM >::Suppress_all_output
protected

If this is set to true then all output from the solver is suppressed. This is protected member data so that the multigrid preconditioner knows whether or not to restore the stream pointer.

Definition at line 642 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::disable_output(), oomph::MGSolver< DIM >::enable_output(), oomph::MGPreconditioner< DIM >::preconditioner_solve(), oomph::MGPreconditioner< DIM >::setup(), and oomph::MGSolver< DIM >::solve().

◆ Suppress_v_cycle_output

template<unsigned DIM>
bool oomph::MGSolver< DIM >::Suppress_v_cycle_output
protected

Indicates whether or not the V-cycle output should be suppressed. Needs to be protected member data for the multigrid preconditioner to know whether or not to output information with each preconditioning step.

Definition at line 637 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::disable_output(), oomph::MGSolver< DIM >::disable_v_cycle_output(), oomph::MGSolver< DIM >::enable_output(), oomph::MGSolver< DIM >::enable_v_cycle_output(), oomph::MGPreconditioner< DIM >::preconditioner_solve(), and oomph::MGSolver< DIM >::solve().

◆ V_cycle_counter

template<unsigned DIM>
unsigned oomph::MGSolver< DIM >::V_cycle_counter
private

Pointer to counter for V-cycles.

Definition at line 727 of file geometric_multigrid.h.

Referenced by oomph::MGSolver< DIM >::iterations(), and oomph::MGSolver< DIM >::solve().

◆ X_mg_vectors_storage

template<unsigned DIM>
Vector<DoubleVector> oomph::MGSolver< DIM >::X_mg_vectors_storage
private

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