general_purpose_block_preconditioners.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 
27 // Include guards
28 #ifndef OOMPH_GENERAL_BLOCK_PRECONDITIONERS
29 #define OOMPH_GENERAL_BLOCK_PRECONDITIONERS
30 
31 
32 // Config header generated by autoconfig
33 #ifdef HAVE_CONFIG_H
34 #include <oomph-lib-config.h>
35 #endif
36 
37 // c++ include
38 // #include<list>
39 
40 // oomph-lib includes
41 #include "matrices.h"
42 // #include "mesh.h"
43 // #include "problem.h"
44 #include "block_preconditioner.h"
45 #include "SuperLU_preconditioner.h"
46 #include "preconditioner_array.h"
47 #include "matrix_vector_product.h"
48 
49 
50 namespace oomph
51 {
52  namespace PreconditionerCreationFunctions
53  {
54  /// Helper function to create a SuperLu preconditioner (for use as
55  /// the default subsididary preconditioner creator in
56  /// GeneralPurposeBlockPreconditioners).
58  {
59  return new SuperLUPreconditioner;
60  }
61  } // namespace PreconditionerCreationFunctions
62 
63 
64  //============================================================================
65  /// Base class for general purpose block preconditioners. Deals with
66  /// setting subsidiary preconditioners and dof to block maps.
67  /// Subsidiary preconditioners can be set in two ways:
68  /// 1) A pointer to a subsidiary preconditioner for block i can be passed
69  /// to set_subsidiary_preconditioner_pt(prec, i).
70  /// 2) A default subsidiary preconditioner can be set up by providing a
71  /// function pointer to a function which creates a preconditioner. During
72  /// setup() all unset subsidiary preconditioner pointers will be filled in
73  /// using this function. By default this uses SuperLU.
74  //============================================================================
75  template<typename MATRIX>
77  {
78  public:
79  /// typedef for a function that allows other preconditioners to be
80  /// employed to solve the subsidiary linear systems.
81  /// The function should return a pointer to the required subsidiary
82  /// preconditioner generated using new. This preconditioner is responsible
83  /// for the destruction of the subsidiary preconditioners.
84  typedef Preconditioner* (*SubsidiaryPreconditionerFctPt)();
85 
86  /// constructor
88  : BlockPreconditioner<MATRIX>(),
90  &PreconditionerCreationFunctions::create_super_lu_preconditioner)
91  {
92  // Make sure that the Gp_mesh_pt container is size zero.
93  Gp_mesh_pt.resize(0);
94  }
95 
96  /// Destructor: clean up memory then delete all subsidiary
97  /// preconditioners.
99  {
100  this->clean_up_memory();
101 
102  for (unsigned j = 0, nj = Subsidiary_preconditioner_pt.size(); j < nj;
103  j++)
104  {
106  }
107  }
108 
109  /// ??ds I think clean_up_memory is supposed to clear out any stuff
110  /// that doesn't need to be stored between solves. Call clean up on any
111  /// non-null subsidiary preconditioners.
112  virtual void clean_up_memory()
113  {
114  // Call clean up in any subsidiary precondtioners that are set.
115  for (unsigned j = 0, nj = Subsidiary_preconditioner_pt.size(); j < nj;
116  j++)
117  {
118  if (Subsidiary_preconditioner_pt[j] != 0)
119  {
120  Subsidiary_preconditioner_pt[j]->clean_up_memory();
121  }
122  }
123 
124  // Clean up the block preconditioner base class stuff
126  }
127 
128  /// Broken copy constructor
130  const GeneralPurposeBlockPreconditioner&) = delete;
131 
132  /// Broken assignment operator
134 
135  /// access function to set the subsidiary preconditioner function.
137  SubsidiaryPreconditionerFctPt sub_prec_fn)
138  {
140  }
141 
142  /// Reset the subsidiary preconditioner function to its default
144  {
147  }
148  /// Set the subsidiary preconditioner to use for block i. The
149  /// subsidiary preconditioner should have been created using new (the
150  /// general purpose block preconditioner will delete it later). If null
151  /// the general purpose block preconditioner will use the
152  /// Subsidiary_preconditioner_creation_function_pt to create the
153  /// preconditioner during setup().
155  const unsigned& i)
156  {
157  // If the vector is currently too small to hold that many
158  // preconditioners then expand it and fill with nulls.
159  if (Subsidiary_preconditioner_pt.size() < i + 1)
160  {
161  Subsidiary_preconditioner_pt.resize(i + 1, 0);
162  }
163  // Note: the size of the vector is checked by
164  // fill_in_subsidiary_preconditioners(..) when we know what size it
165  // should be.
166 
167  // I'm assuming that the number of preconditioners is always "small"
168  // compared to Jacobian size, so a resize doesn't waste much time.
169 
170  // Put the pointer in the vector
172  }
173 
174  /// Get the subsidiary precondtioner pointer in block i (is
175  /// allowed to be null if not yet set).
177  {
179  }
180 
181  /// Specify a DOF to block map
182  void set_dof_to_block_map(Vector<unsigned>& dof_to_block_map)
183  {
184  Dof_to_block_map = dof_to_block_map;
185  }
186 
187  /// Adds a mesh to be used by the
188  /// block preconditioning framework for classifying DOF types. Optional
189  /// boolean argument (default: false) allows the mesh to contain multiple
190  /// element types.
191  void add_mesh(const Mesh* mesh_pt,
192  const bool& allow_multiple_element_type_in_mesh = false)
193  {
194 #ifdef PARANOID
195  // Check that the mesh pointer is not null.
196  if (mesh_pt == 0)
197  {
198  std::ostringstream err_msg;
199  err_msg << "The mesh_pt is null, please point it to a mesh.\n";
200  throw OomphLibError(
201  err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
202  }
203 #endif
204  // Push back the mesh pointer and the boolean in a pair.
205  Gp_mesh_pt.push_back(
206  std::make_pair(mesh_pt, allow_multiple_element_type_in_mesh));
207  }
208 
209  /// Returns the number of meshes currently set in the
210  /// GeneralPurposeBlockPreconditioner base class.
211  unsigned gp_nmesh()
212  {
213  return Gp_mesh_pt.size();
214  }
215 
216  protected:
217  /// Set the mesh in the block preconditioning framework.
219  {
220  const unsigned nmesh = gp_nmesh();
221 #ifdef PARANOID
222  if (nmesh == 0)
223  {
224  std::ostringstream err_msg;
225  err_msg << "There are no meshes set.\n"
226  << "Have you remembered to call add_mesh(...)?\n";
227  throw OomphLibError(
228  err_msg.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
229  }
230 #endif
231 
232  this->set_nmesh(nmesh);
233  for (unsigned mesh_i = 0; mesh_i < nmesh; mesh_i++)
234  {
235  this->set_mesh(
236  mesh_i, Gp_mesh_pt[mesh_i].first, Gp_mesh_pt[mesh_i].second);
237  }
238  }
239 
240  /// Modified block setup for general purpose block preconditioners
242  {
243  if (Dof_to_block_map.size() > 0)
244  {
246  }
247  else
248  {
250  }
251  }
252 
253  /// Create any subsidiary preconditioners needed. Usually
254  /// nprec_needed = nblock_types, except for the ExactBlockPreconditioner
255  /// which only requires one preconditioner.
256  void fill_in_subsidiary_preconditioners(const unsigned& nprec_needed)
257  {
258  // If it's empty then fill it in with null pointers.
259  if (Subsidiary_preconditioner_pt.empty())
260  {
261  Subsidiary_preconditioner_pt.assign(nprec_needed, 0);
262  }
263  else
264  {
265  // Otherwise check we have the right number of them
266 #ifdef PARANOID
267  if (Subsidiary_preconditioner_pt.size() != nprec_needed)
268  {
269  using namespace StringConversion;
270  std::string error_msg = "Wrong number of precondtioners in";
271  error_msg += "Subsidiary_preconditioner_pt, should have ";
272  error_msg += to_string(nprec_needed) + " but we actually have ";
273  error_msg += to_string(Subsidiary_preconditioner_pt.size());
274  throw OomphLibError(
275  error_msg, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
276  }
277 #endif
278  }
279 
280 
281  // Now replace any null pointers with new preconditioners
282  for (unsigned j = 0, nj = Subsidiary_preconditioner_pt.size(); j < nj;
283  j++)
284  {
285  if (Subsidiary_preconditioner_pt[j] == 0)
286  {
288  (*Subsidiary_preconditioner_creation_function_pt)();
289  }
290  }
291  }
292 
293  /// List of preconditioners to use for the blocks to be solved.
295 
296  /// Function to create any subsidiary preconditioners not set in
297  /// Subsidiary_preconditioner_pt.
300 
301  private:
302  /// the set of dof to block maps for this preconditioner
304 
305  /// Vector of mesh pointers and a boolean indicating if we allow multiple
306  /// element types in the same mesh.
308  };
309 
310 
311  //=============================================================================
312  /// Block diagonal preconditioner. By default SuperLU is used to solve
313  /// the subsidiary systems, but other preconditioners can be used by setting
314  /// them using passing a pointer to a function of type
315  /// SubsidiaryPreconditionerFctPt to the method
316  /// subsidiary_preconditioner_function_pt().
317  //=============================================================================
318  template<typename MATRIX>
320  : public GeneralPurposeBlockPreconditioner<MATRIX>
321  {
322  public:
323  /// constructor - when the preconditioner is used as a master preconditioner
325  {
326  // by default we do not use two level parallelism
327  Use_two_level_parallelisation = false;
328 
329  // null the Preconditioner array pt
330  Preconditioner_array_pt = 0;
331 
332  // Don't doc by default
333  Doc_time_during_preconditioner_solve = false;
334  }
335 
336  /// Destructor - delete the preconditioner matrices
338  {
339  this->clean_up_memory();
340  }
341 
342  /// clean up the memory
343  virtual void clean_up_memory()
344  {
345  if (Use_two_level_parallelisation)
346  {
347  delete Preconditioner_array_pt;
348  Preconditioner_array_pt = 0;
349  }
350 
351  // Clean up the base class too
353  }
354 
355  /// Broken copy constructor
357 
358  /// Broken assignment operator
359  void operator=(const BlockDiagonalPreconditioner&) = delete;
360 
361  /// Apply preconditioner to r
363 
364  /// Setup the preconditioner
365  virtual void setup();
366 
367  /// Use two level parallelisation
369  {
370 #ifndef OOMPH_HAS_MPI
371  throw OomphLibError("Cannot do any parallelism since we don't have MPI.",
372  OOMPH_CURRENT_FUNCTION,
373  OOMPH_EXCEPTION_LOCATION);
374 #endif
375  Use_two_level_parallelisation = true;
376  }
377 
378  /// Don't use two-level parallelisation
380  {
381  Use_two_level_parallelisation = false;
382  }
383 
384  /// Enable Doc timings in application of block sub-preconditioners
386  {
387  Doc_time_during_preconditioner_solve = true;
388  }
389 
390  /// Disable Doc timings in application of block sub-preconditioners
392  {
393  Doc_time_during_preconditioner_solve = false;
394  }
395 
396  void fill_in_subsidiary_preconditioners(const unsigned& nprec_needed)
397  {
398 #ifdef PARANOID
399  if ((Use_two_level_parallelisation) &&
400  !this->Subsidiary_preconditioner_pt.empty())
401  {
402  std::string err_msg =
403  "Two level parallelism diagonal block preconditioners cannot have";
404  err_msg +=
405  " any preset preconditioners (due to weird memory management";
406  err_msg += " in the PreconditionerArray, you could try fixing it).";
407  throw OomphLibError(
408  err_msg, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
409  }
410 #endif
411 
412  // Now call the real function
414  MATRIX>::fill_in_subsidiary_preconditioners(nprec_needed);
415  }
416 
417  protected:
418  /// This is a helper function to allow us to implement AntiDiagonal
419  /// preconditioner by only changing this function. Get the second index
420  /// for block number i. Obviously for a diagonal preconditioner we want
421  /// the blocks (i,i), (for anti diagonal we will want blocks (i, nblock -
422  /// i), see that class).
423  virtual unsigned get_other_diag_ds(const unsigned& i,
424  const unsigned& nblock) const
425  {
426  return i;
427  }
428 
429 
430  private:
431  /// pointer for the PreconditionerArray
433 
434  /// Use two level parallelism using the PreconditionerArray
436 
437  /// Doc timings in application of block sub-preconditioners?
439  };
440 
441 
442  /// ////////////////////////////////////////////////////////////////////////////
443  /// ////////////////////////////////////////////////////////////////////////////
444  /// ////////////////////////////////////////////////////////////////////////////
445 
446 
447  //=============================================================================
448  /// General purpose block triangular preconditioner
449  /// By default this is Upper triangular.
450  /// By default SuperLUPreconditioner (or SuperLUDistPreconditioner) is used to
451  /// solve the subsidiary systems, but other preconditioners can be used by
452  /// setting them using passing a pointer to a function of type
453  /// SubsidiaryPreconditionerFctPt to the method
454  /// subsidiary_preconditioner_function_pt().
455  //=============================================================================
456  template<typename MATRIX>
458  : public GeneralPurposeBlockPreconditioner<MATRIX>
459  {
460  public:
461  /// Constructor. (By default this preconditioner is upper triangular).
464  {
465  // default to upper triangular
466  Upper_triangular = true;
467  }
468 
469  /// Destructor - delete the preconditioner matrices
471  {
472  this->clean_up_memory();
473  }
474 
475  /// clean up the memory
476  virtual void clean_up_memory()
477  {
478  // Delete anything in Off_diagonal_matrix_vector_products
479  for (unsigned i = 0, ni = Off_diagonal_matrix_vector_products.nrow();
480  i < ni;
481  i++)
482  {
483  for (unsigned j = 0, nj = Off_diagonal_matrix_vector_products.ncol();
484  j < nj;
485  j++)
486  {
487  delete Off_diagonal_matrix_vector_products(i, j);
488  Off_diagonal_matrix_vector_products(i, j) = 0;
489  }
490  }
491 
492  // Clean up the base class too
494  }
495 
496  /// Broken copy constructor
498  delete;
499 
500  /// Broken assignment operator
502 
503  /// Apply preconditioner to r
505 
506  /// Setup the preconditioner
507  void setup();
508 
509  /// Use as an upper triangular preconditioner
511  {
512  Upper_triangular = true;
513  }
514 
515  /// Use as a lower triangular preconditioner
517  {
518  Upper_triangular = false;
519  }
520 
521  private:
522  /// Matrix of matrix vector product operators for the off diagonals
524 
525  /// Boolean indicating upper or lower triangular
527  };
528 
529 
530  /// ////////////////////////////////////////////////////////////////////////////
531  /// ////////////////////////////////////////////////////////////////////////////
532  /// ////////////////////////////////////////////////////////////////////////////
533 
534 
535  //=============================================================================
536  /// Exact block preconditioner - block preconditioner assembled from all
537  /// blocks associated with the preconditioner and solved by SuperLU.
538  //=============================================================================
539  template<typename MATRIX>
541  : public GeneralPurposeBlockPreconditioner<MATRIX>
542  {
543  public:
544  /// constructor
546 
547  /// Destructor
549 
550  /// Broken copy constructor
552 
553  /// Broken assignment operator
554  void operator=(const ExactBlockPreconditioner&) = delete;
555 
556  /// Apply preconditioner to r
558 
559  /// Setup the preconditioner
560  void setup();
561 
562  /// Access for the preconditioner pointer used to solve the
563  /// system (stored in the vector of pointers in the base class);
565  {
566  return this->Subsidiary_preconditioner_pt[0];
567  }
568  };
569 
570 
571  // =================================================================
572  /// Block "anti-diagonal" preconditioner, i.e. same as block
573  /// diagonal but along the other diagonal of the matrix (top-right to
574  /// bottom-left).
575  // =================================================================
576  template<typename MATRIX>
578  : public BlockDiagonalPreconditioner<MATRIX>
579  {
580  protected:
581  /// This is a helper function to allow us to implement BlockAntiDiagonal
582  /// using BlockDiagonal preconditioner and only changing this
583  /// function. Get the second index for block number i. Obviously for a
584  /// diagonal preconditioner we want the blocks (i,i). For anti diagonal
585  /// we will want blocks (i, nblock - i - 1).
586  unsigned get_other_diag_ds(const unsigned& i, const unsigned& nblock) const
587  {
588  return nblock - i - 1;
589  }
590  };
591 
592 
593  // =================================================================
594  /// Preconditioner that doesn't actually do any preconditioning, it just
595  /// allows access to the Jacobian blocks. This is pretty hacky but oh well..
596  // =================================================================
597  template<typename MATRIX>
599  : public GeneralPurposeBlockPreconditioner<MATRIX>
600  {
601  public:
602  /// Constructor
604 
605  /// Destructor
607 
608  /// Broken copy constructor
610 
611  /// Broken assignment operator
612  void operator=(const DummyBlockPreconditioner&) = delete;
613 
614  /// Apply preconditioner to r (just copy r to z).
616  {
617  z.build(r);
618  }
619 
620  /// Setup the preconditioner
621  void setup()
622  {
623  // Set up the block look up schemes
624  this->block_setup();
625  }
626  };
627 
628 } // namespace oomph
629 #endif
cstr elem_len * i
Definition: cfortran.h:603
Block "anti-diagonal" preconditioner, i.e. same as block diagonal but along the other diagonal of the...
unsigned get_other_diag_ds(const unsigned &i, const unsigned &nblock) const
This is a helper function to allow us to implement BlockAntiDiagonal using BlockDiagonal precondition...
Block diagonal preconditioner. By default SuperLU is used to solve the subsidiary systems,...
void disable_doc_time_during_preconditioner_solve()
Disable Doc timings in application of block sub-preconditioners.
virtual ~BlockDiagonalPreconditioner()
Destructor - delete the preconditioner matrices.
void fill_in_subsidiary_preconditioners(const unsigned &nprec_needed)
void enable_two_level_parallelisation()
Use two level parallelisation.
bool Doc_time_during_preconditioner_solve
Doc timings in application of block sub-preconditioners?
void operator=(const BlockDiagonalPreconditioner &)=delete
Broken assignment operator.
void enable_doc_time_during_preconditioner_solve()
Enable Doc timings in application of block sub-preconditioners.
void disable_two_level_parallelisation()
Don't use two-level parallelisation.
BlockDiagonalPreconditioner(const BlockDiagonalPreconditioner &)=delete
Broken copy constructor.
BlockDiagonalPreconditioner()
constructor - when the preconditioner is used as a master preconditioner
virtual unsigned get_other_diag_ds(const unsigned &i, const unsigned &nblock) const
This is a helper function to allow us to implement AntiDiagonal preconditioner by only changing this ...
PreconditionerArray * Preconditioner_array_pt
pointer for the PreconditionerArray
bool Use_two_level_parallelisation
Use two level parallelism using the PreconditionerArray.
Block Preconditioner base class. The block structure of the overall problem is determined from the Me...
unsigned nmesh() const
Return the number of meshes in Mesh_pt.
const Mesh * mesh_pt(const unsigned &i) const
Access to i-th mesh (of the various meshes that contain block preconditionable elements of the same n...
void clear_block_preconditioner_base()
Clears all BlockPreconditioner data. Called by the destructor and the block_setup(....
void set_nmesh(const unsigned &n)
Specify the number of meshes required by this block preconditioner. Note: elements in different meshe...
virtual void block_setup()
Determine the size of the matrix blocks and setup the lookup schemes relating the global degrees of f...
void set_mesh(const unsigned &i, const Mesh *const mesh_pt, const bool &allow_multiple_element_type_in_mesh=false)
Set the i-th mesh for this block preconditioner. Note: The method set_nmesh(...) must be called befor...
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
virtual ~BlockTriangularPreconditioner()
Destructor - delete the preconditioner matrices.
DenseMatrix< MatrixVectorProduct * > Off_diagonal_matrix_vector_products
Matrix of matrix vector product operators for the off diagonals.
void upper_triangular()
Use as an upper triangular preconditioner.
BlockTriangularPreconditioner()
Constructor. (By default this preconditioner is upper triangular).
void operator=(const BlockTriangularPreconditioner &)=delete
Broken assignment operator.
BlockTriangularPreconditioner(const BlockTriangularPreconditioner &)=delete
Broken copy constructor.
void lower_triangular()
Use as a lower triangular preconditioner.
bool Upper_triangular
Boolean indicating upper or lower triangular.
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
Definition: matrices.h:386
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
void build(const DoubleVector &old_vector)
Just copys the argument DoubleVector.
Preconditioner that doesn't actually do any preconditioning, it just allows access to the Jacobian bl...
DummyBlockPreconditioner(const DummyBlockPreconditioner &)=delete
Broken copy constructor.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r (just copy r to z).
void operator=(const DummyBlockPreconditioner &)=delete
Broken assignment operator.
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
void operator=(const ExactBlockPreconditioner &)=delete
Broken assignment operator.
ExactBlockPreconditioner(const ExactBlockPreconditioner &)=delete
Broken copy constructor.
Preconditioner *& preconditioner_pt()
Access for the preconditioner pointer used to solve the system (stored in the vector of pointers in t...
Base class for general purpose block preconditioners. Deals with setting subsidiary preconditioners a...
void operator=(const GeneralPurposeBlockPreconditioner &)=delete
Broken assignment operator.
Vector< unsigned > Dof_to_block_map
the set of dof to block maps for this preconditioner
void set_dof_to_block_map(Vector< unsigned > &dof_to_block_map)
Specify a DOF to block map.
virtual void clean_up_memory()
??ds I think clean_up_memory is supposed to clear out any stuff that doesn't need to be stored betwee...
virtual ~GeneralPurposeBlockPreconditioner()
Destructor: clean up memory then delete all subsidiary preconditioners.
Preconditioner *(* SubsidiaryPreconditionerFctPt)()
typedef for a function that allows other preconditioners to be employed to solve the subsidiary linea...
SubsidiaryPreconditionerFctPt Subsidiary_preconditioner_creation_function_pt
Function to create any subsidiary preconditioners not set in Subsidiary_preconditioner_pt.
Vector< Preconditioner * > Subsidiary_preconditioner_pt
List of preconditioners to use for the blocks to be solved.
Vector< std::pair< const Mesh *, bool > > Gp_mesh_pt
Vector of mesh pointers and a boolean indicating if we allow multiple element types in the same mesh.
Preconditioner * subsidiary_preconditioner_pt(const unsigned &i) const
Get the subsidiary precondtioner pointer in block i (is allowed to be null if not yet set).
void reset_subsidiary_preconditioner_function_to_default()
Reset the subsidiary preconditioner function to its default.
void add_mesh(const Mesh *mesh_pt, const bool &allow_multiple_element_type_in_mesh=false)
Adds a mesh to be used by the block preconditioning framework for classifying DOF types....
void gp_preconditioner_set_all_meshes()
Set the mesh in the block preconditioning framework.
GeneralPurposeBlockPreconditioner(const GeneralPurposeBlockPreconditioner &)=delete
Broken copy constructor.
void fill_in_subsidiary_preconditioners(const unsigned &nprec_needed)
Create any subsidiary preconditioners needed. Usually nprec_needed = nblock_types,...
void set_subsidiary_preconditioner_pt(Preconditioner *prec, const unsigned &i)
Set the subsidiary preconditioner to use for block i. The subsidiary preconditioner should have been ...
void gp_preconditioner_block_setup()
Modified block setup for general purpose block preconditioners.
void set_subsidiary_preconditioner_function(SubsidiaryPreconditionerFctPt sub_prec_fn)
access function to set the subsidiary preconditioner function.
unsigned gp_nmesh()
Returns the number of meshes currently set in the GeneralPurposeBlockPreconditioner base class.
A general mesh class.
Definition: mesh.h:67
An OomphLibError object which should be thrown when an run-time error is encountered....
PreconditionerArray - NOTE - first implementation, a number of assumptions / simplifications were mad...
Preconditioner base class. Gives an interface to call all other preconditioners through and stores th...
virtual void preconditioner_solve(const DoubleVector &r, DoubleVector &z)=0
Apply the preconditioner. Pure virtual generic interface function. This method should apply the preco...
virtual void setup()=0
Setup the preconditioner. Pure virtual generic interface function.
An interface to allow SuperLU to be used as an (exact) Preconditioner.
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
Preconditioner * create_super_lu_preconditioner()
Helper function to create a SuperLu preconditioner (for use as the default subsididary preconditioner...
std::string to_string(T object, unsigned float_precision=8)
Conversion function that should work for anything with operator<< defined (at least all basic types).
//////////////////////////////////////////////////////////////////// ////////////////////////////////...