biharmonic_preconditioner.h
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2023 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Include guards
27 #ifndef OOMPH_BIHARMONIC_PRECONDITIONER_HEADER
28 #define OOMPH_BIHARMONIC_PRECONDITIONER_HEADER
29 
30 
31 // Config header generated by autoconfig
32 #ifdef HAVE_CONFIG_H
33 #include <oomph-lib-config.h>
34 #endif
35 
36 #include "../generic/preconditioner.h"
37 #include "../generic/block_preconditioner.h"
38 #include "../generic/hijacked_elements.h"
39 #include "biharmonic_elements.h"
40 #include "../meshes/hermite_element_quad_mesh.template.h"
41 #include "../generic/SuperLU_preconditioner.h"
42 #include "../generic/general_purpose_preconditioners.h"
43 
44 #ifdef OOMPH_HAS_HYPRE
45 #include "../generic/hypre_solver.h"
46 #endif
47 
48 namespace oomph
49 {
50 #ifdef OOMPH_HAS_HYPRE
51 
52  //=============================================================================
53  // defaults settings for the Hypre solver (AMG) when used as the approximate
54  // linear solver for the Schur complement (non-compound) linear subsidiary
55  // linear systems
56  //=============================================================================
57  namespace Biharmonic_schur_complement_Hypre_defaults
58  {
59  /// smoother type - Gauss Seidel: 1
60  extern unsigned AMG_smoother;
61 
62  /// amg coarsening strategy: classical Ruge Stueben: 1
63  extern unsigned AMG_coarsening;
64 
65  /// number of V cycles: 2
66  extern unsigned N_cycle;
67 
68  /// amg strength parameter: 0.25 -- optimal for 2d
69  extern double AMG_strength;
70 
71  /// jacobi damping -- hierher not used 0.1;
72  extern double AMG_jacobi_damping;
73 
74  /// Amg smoother iterations: 2
75  extern unsigned AMG_smoother_iterations;
76 
77  /// set the defaults
78  extern void set_defaults(HyprePreconditioner* hypre_prec_pt);
79 
80  } // namespace Biharmonic_schur_complement_Hypre_defaults
81 #endif
82 
83 
84  //=============================================================================
85  /// Biharmonic Preconditioner - for two dimensional problems
86  //=============================================================================
87  class BiharmonicPreconditioner : public BlockPreconditioner<CRDoubleMatrix>
88  {
89  public:
90  /// Constructor - by default inexact preconditioning is used
92  {
93  // initialise both preconditioners to zero
97 
98  // by default we use the inexact biharmonic preconditioner
99  // and if possible use Hypre preconditioner
100 #ifdef OOMPH_HAS_HYPRE
102 #else
104 #endif
105 
106  // size mesh pt correctly
107  this->set_nmesh(1);
109  }
110 
111  /// destructor - cleans up preconditioners and delete matrices
113  {
114  this->clean_up_memory();
115  }
116 
117  // delete the subsidiary preconditioners and memory
119  {
120  // delete the sub preconditioners
124  }
125 
126  /// Broken copy constructor
128 
129  /// Broken assignment operator
130  void operator=(const BiharmonicPreconditioner&) = delete;
131 
132  /// Setup the preconditioner
133  void setup();
134 
135  /// Apply preconditioner to r
137 
138  /// Access function to the preconditioner type \n
139  /// + 0 : exact BBD \n
140  /// + 1 : inexact BBD w/ SuperLU \n
141  /// + 2 : inexact BBD w/ AMG (Hypre Boomer AMG)
142  /// + 3 : block diagonal (3x3)+(1x1)
144  {
145  return Preconditioner_type;
146  }
147 
148  /// Access function to the mesh pt for the bulk elements. The mesh
149  /// should only contain BiharmonicElement<2> and
150  /// Hijacked<BiharmonicElement<2> > elements
152  {
153  return Bulk_element_mesh_pt;
154  }
155 
156  private:
157  /// preconditioner type \n
158  /// + 0 : exact BBD \n
159  /// + 1 : inexact BBD w/ SuperLU \n
160  /// + 2 : inexact BBD w/ AMG
161  /// + 3 : block diagonal (3x3)+(1x1)
163 
164  /// Exact Preconditioner Pointer
166 
167  /// Inexact Preconditioner Pointer
169 
170  /// Preconditioner the diagonal block associated with hijacked elements
172 
173  /// the bulk element mesh pt
175  };
176 
177 
178  //=============================================================================
179  /// Sub Biharmonic Preconditioner - an exact preconditioner for the
180  /// 3x3 top left hand corner sub block matrix.
181  /// Used as part of the BiharmonicPreconditioner<MATRIX> .
182  /// By default this uses the BBD (block-bordered-diagonal/arrow-shaped)
183  /// preconditioner; can also switch to full BD version (in which case
184  /// all the 3x3 blocks are retained)
185  //=============================================================================
187  : public BlockPreconditioner<CRDoubleMatrix>
188  {
189  public:
190  /// Constructor - for a preconditioner acting as a sub preconditioner
192  const bool& retain_all_blocks = false)
193  : Retain_all_blocks(retain_all_blocks)
194  {
195  // Block mapping for ExactSubBiharmonicPreconditioner
196  Vector<unsigned> block_lookup(3);
197  block_lookup[0] = 0;
198  block_lookup[1] = 1;
199  block_lookup[2] = 2;
200 
201  // set as subsidiary block preconditioner
202  this->turn_into_subsidiary_block_preconditioner(master_prec_pt,
203  block_lookup);
204 
205  // null the Sub preconditioner pt
207  }
208 
209  /// destructor deletes the exact preconditioner
211  {
212  this->clean_up_memory();
213  }
214 
215  /// delete the subsidiary preconditioner pointer
217  {
218  delete Sub_preconditioner_pt;
220  }
221 
222  /// Broken copy constructor
224  delete;
225 
226  /// Broken assignment operator
228 
229  /// Setup the preconditioner
230  void setup();
231 
232  /// Apply preconditioner to r
234 
235  // private:
236 
237  // Pointer to the sub preconditioner
239 
240 
241  /// Boolean indicating that all blocks are to be retained (defaults to
242  /// false)
244  };
245 
246  //=============================================================================
247  /// SubBiharmonic Preconditioner - an inexact preconditioner for the
248  /// 3x3 top left hand corner sub block matrix.
249  /// Used as part of the BiharmonicPreconditioner<MATRIX>
250  //=============================================================================
252  : public BlockPreconditioner<CRDoubleMatrix>
253  {
254  public:
255  /// Constructor for the inexact block preconditioner. \n
256  /// This a helper class for BiharmonicPreconditioner and cannot be used
257  /// as a standalone preconditioner. \n
258  /// master_prec_pt is the pointer to the master BiharmonicPreconditioner.
260  const bool use_amg)
261  {
262  // Block mapping for ExactSubBiharmonicPreconditioner
263  Vector<unsigned> block_lookup(3);
264  block_lookup[0] = 0;
265  block_lookup[1] = 1;
266  block_lookup[2] = 2;
267 
268  // set as subsidiary block preconditioner
269  this->turn_into_subsidiary_block_preconditioner(master_prec_pt,
270  block_lookup);
271 
272  // zero all the preconditioner pt
276 
277  // store the preconditioner type
278  Use_amg = use_amg;
279  }
280 
281  /// destructor - just calls this->clean_up_memory()
283  {
284  this->clean_up_memory();
285  }
286 
287  /// clean the memory
289  {
290  // delete the sub preconditioner pt
291  delete S_00_preconditioner_pt;
293 
294  // delete the lumped preconditioners
299 
300  // Number of block types
301  unsigned n = Matrix_of_block_pointers.nrow();
302 
303  // delete the block matrices
304  for (unsigned i = 0; i < n; i++)
305  {
306  for (unsigned j = 0; j < n; j++)
307  {
308  if (Matrix_of_block_pointers(i, j) != 0)
309  {
310  delete Matrix_of_block_pointers(i, j);
311  Matrix_of_block_pointers(i, j) = 0;
312  }
313  }
314  }
315  }
316 
317  /// Broken copy constructor
319  const InexactSubBiharmonicPreconditioner&) = delete;
320 
321  /// Broken assignment operator
323 
324  /// Setup the preconditioner
325  void setup();
326 
327  /// Apply preconditioner to r
329 
330  // private:
331 
332  /// Computes the inexact schur complement of the block J_00 using
333  /// lumping as an approximate inverse of block J_11 and J_22
334  /// (where possible)
336 
337  /// Pointer to the S00 Schur Complement preconditioner
339 
340  /// Preconditioner for storing the lumped J_11 matrix
343 
344  /// Preconditioner for storing the lumped J_22 matrix
347 
348  ///
350 
351  ///
353 
354  /// booean indicating whether (Hypre BoomerAMG) AMG should be used
355  /// to solve the S00 subsidiary linear system.
356  unsigned Use_amg;
357  };
358 } // namespace oomph
359 #endif
cstr elem_len * i
Definition: cfortran.h:603
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.
Definition: matrices.h:888
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
Definition: matrices.h:386
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
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.
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
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...
~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 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.
A general mesh class.
Definition: mesh.h:67
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
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
//////////////////////////////////////////////////////////////////// ////////////////////////////////...