27 #ifndef OOMPH_BIHARMONIC_PRECONDITIONER_HEADER
28 #define OOMPH_BIHARMONIC_PRECONDITIONER_HEADER
33 #include <oomph-lib-config.h>
36 #include "../generic/preconditioner.h"
37 #include "../generic/block_preconditioner.h"
38 #include "../generic/hijacked_elements.h"
40 #include "../meshes/hermite_element_quad_mesh.template.h"
41 #include "../generic/SuperLU_preconditioner.h"
42 #include "../generic/general_purpose_preconditioners.h"
44 #ifdef OOMPH_HAS_HYPRE
45 #include "../generic/hypre_solver.h"
50 #ifdef OOMPH_HAS_HYPRE
57 namespace Biharmonic_schur_complement_Hypre_defaults
78 extern void set_defaults(HyprePreconditioner* hypre_prec_pt);
100 #ifdef OOMPH_HAS_HYPRE
192 const bool& retain_all_blocks =
false)
304 for (
unsigned i = 0;
i < n;
i++)
306 for (
unsigned j = 0; j < n; j++)
Biharmonic Preconditioner - for two dimensional problems.
BiharmonicPreconditioner()
Constructor - by default inexact preconditioning is used.
Preconditioner * Hijacked_sub_block_preconditioner_pt
Preconditioner the diagonal block associated with hijacked elements.
Preconditioner * Sub_preconditioner_2_pt
Inexact Preconditioner Pointer.
void clean_up_memory()
Clean up memory (empty). Generic interface function.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r.
BiharmonicPreconditioner(const BiharmonicPreconditioner &)=delete
Broken copy constructor.
~BiharmonicPreconditioner()
destructor - cleans up preconditioners and delete matrices
void operator=(const BiharmonicPreconditioner &)=delete
Broken assignment operator.
unsigned & preconditioner_type()
Access function to the preconditioner type .
void setup()
Setup the preconditioner.
unsigned Preconditioner_type
preconditioner type
Mesh * Bulk_element_mesh_pt
the bulk element mesh pt
Mesh *& bulk_element_mesh_pt()
Access function to the mesh pt for the bulk elements. The mesh should only contain BiharmonicElement<...
Preconditioner * Sub_preconditioner_1_pt
Exact Preconditioner Pointer.
Block Preconditioner base class. The block structure of the overall problem is determined from the Me...
void turn_into_subsidiary_block_preconditioner(BlockPreconditioner< CRDoubleMatrix > *master_block_prec_pt, const Vector< unsigned > &doftype_in_master_preconditioner_coarse)
Function to turn this preconditioner into a subsidiary preconditioner that operates within a bigger "...
void set_nmesh(const unsigned &n)
Specify the number of meshes required by this block preconditioner. Note: elements in different meshe...
A class for compressed row matrices. This is a distributable object.
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
A vector in the mathematical sense, initially developed for linear algebra type applications....
Sub Biharmonic Preconditioner - an exact preconditioner for the 3x3 top left hand corner sub block ma...
ExactSubBiharmonicPreconditioner(BiharmonicPreconditioner *master_prec_pt, const bool &retain_all_blocks=false)
Constructor - for a preconditioner acting as a sub preconditioner.
ExactSubBiharmonicPreconditioner(const ExactSubBiharmonicPreconditioner &)=delete
Broken copy constructor.
void clean_up_memory()
delete the subsidiary preconditioner pointer
void operator=(const ExactSubBiharmonicPreconditioner &)=delete
Broken assignment operator.
Preconditioner * Sub_preconditioner_pt
bool Retain_all_blocks
Boolean indicating that all blocks are to be retained (defaults to false)
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r.
~ExactSubBiharmonicPreconditioner()
destructor deletes the exact preconditioner
void setup()
Setup the preconditioner.
SubBiharmonic Preconditioner - an inexact preconditioner for the 3x3 top left hand corner sub block m...
DenseMatrix< CRDoubleMatrix * > Matrix_of_block_pointers
void compute_inexact_schur_complement()
Computes the inexact schur complement of the block J_00 using lumping as an approximate inverse of bl...
Preconditioner * S_00_preconditioner_pt
Pointer to the S00 Schur Complement preconditioner.
unsigned Use_amg
booean indicating whether (Hypre BoomerAMG) AMG should be used to solve the S00 subsidiary linear sys...
InexactSubBiharmonicPreconditioner(BiharmonicPreconditioner *master_prec_pt, const bool use_amg)
Constructor for the inexact block preconditioner. This a helper class for BiharmonicPreconditioner a...
void clean_up_memory()
clean the memory
~InexactSubBiharmonicPreconditioner()
destructor - just calls this->clean_up_memory()
void operator=(const InexactSubBiharmonicPreconditioner &)=delete
Broken assignment operator.
MatrixBasedLumpedPreconditioner< CRDoubleMatrix > * Lumped_J_22_preconditioner_pt
Preconditioner for storing the lumped J_22 matrix.
InexactSubBiharmonicPreconditioner(const InexactSubBiharmonicPreconditioner &)=delete
Broken copy constructor.
void setup()
Setup the preconditioner.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r.
MatrixBasedLumpedPreconditioner< CRDoubleMatrix > * Lumped_J_11_preconditioner_pt
Preconditioner for storing the lumped J_11 matrix.
Matrix-based lumped preconditioner.
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
unsigned N_cycle
number of V cycles: 2
double AMG_jacobi_damping
jacobi damping – hierher not used 0.1
unsigned AMG_smoother
smoother type - Gauss Seidel: 1
void set_defaults(HyprePreconditioner *hypre_prec_pt)
set the defaults
unsigned AMG_coarsening
amg coarsening strategy: classical Ruge Stueben: 1
double AMG_strength
amg strength parameter: 0.25 – optimal for 2d
unsigned AMG_smoother_iterations
amg smoother iterations
//////////////////////////////////////////////////////////////////// ////////////////////////////////...