biharmonic_preconditioner.cc
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 // Config header generated by autoconfig
27 #ifdef HAVE_CONFIG_H
28 #include <oomph-lib-config.h>
29 #endif
30 
31 
32 // oomph-lib includes
34 
35 
36 namespace oomph
37 {
38 #ifdef OOMPH_HAS_HYPRE
39 
40  //=============================================================================
41  // defaults settings for the Hypre solver (AMG) when used as the approximate
42  // linear solver for the Schur complement (non-compound) linear subsidiary
43  // linear systems
44  //=============================================================================
45  namespace Biharmonic_schur_complement_Hypre_defaults
46  {
47  /// smoother type - Gauss Seidel: 1
48  unsigned AMG_smoother = 1;
49 
50  /// amg coarsening strategy: classical Ruge Stueben: 1
51  unsigned AMG_coarsening = 1;
52 
53  /// number of V cycles: 2
54  unsigned N_cycle = 2;
55 
56  /// amg strength parameter: 0.25 -- optimal for 2d
57  double AMG_strength = 0.25;
58 
59  /// jacobi damping -- hierher not used 0.1
60  double AMG_jacobi_damping = 0.1;
61 
62  /// amg smoother iterations
64 
65  /// set the defaults
66  void set_defaults(HyprePreconditioner* hypre_prec_pt)
67  {
68  // use AMG preconditioner
69  hypre_prec_pt->hypre_method() = HypreSolver::BoomerAMG;
70 
71  // Smoother types
72  hypre_prec_pt->amg_simple_smoother() = AMG_smoother;
73 
74  // jacobi damping
75  // hypre_prec_pt->amg_damping() = AMG_jacobi_damping;
76 
77  // coarsening stategy
78  hypre_prec_pt->amg_coarsening() = AMG_coarsening;
79 
80  oomph_info << "Current number of v cycles: "
81  << hypre_prec_pt->amg_iterations() << std::endl;
82 
83  // number of v-cycles
84  hypre_prec_pt->amg_iterations() = N_cycle;
85 
86  oomph_info << "Re-assigned number of v cycles: "
87  << hypre_prec_pt->amg_iterations() << std::endl;
88 
89  // strength parameter
90  hypre_prec_pt->amg_strength() = AMG_strength;
91 
92  // hierher new
93  oomph_info << "Current number of amg smoother iterations: "
94  << hypre_prec_pt->amg_smoother_iterations() << std::endl;
95 
97 
98  oomph_info << "Re-assigned number of amg smoother iterations: "
99  << hypre_prec_pt->amg_smoother_iterations() << std::endl;
100  }
101  } // namespace Biharmonic_schur_complement_Hypre_defaults
102 #endif
103 
104  //===========================================================================
105  /// setup for the biharmonic preconditioner
106  //===========================================================================
108  {
109  // clean up
110  this->clean_up_memory();
111 
112  // paranoid check that teh bulk element mesh has been set
113 #ifdef PARANOID
114  if (Bulk_element_mesh_pt == 0)
115  {
116  std::ostringstream error_message;
117  error_message << "The bulk element mesh has not been passed to "
118  "bulk_element_mesh_pt()";
119  throw OomphLibError(
120  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
121  }
122 #endif
123 
124  // setup the mesh
125  this->set_mesh(0, Bulk_element_mesh_pt);
126 
127  // setup the blocks look up schemes
128  this->block_setup();
129 
130  // determine whether this preconditioner has 4 or 5 block types and set
131  // Nblock_types if neccessary
132  // unsigned n_row = this->master_nrow();
133  // bool nblock_type_check = true;
134  // for (unsigned i = 0; i < n_row; i++)
135  // {
136  // if (this->block_number(i) == 4) { nblock_type_check = false; }
137  // }
138  // if (nblock_type_check) { Nblock_types = 4; }
139  //
140 
141  // check the preconditioner type is acceptable
142 #ifdef PARANOID
143  if (Preconditioner_type != 0 && Preconditioner_type != 1 &&
145  {
146  std::ostringstream error_message;
147  error_message << "Preconditioner_type must be equal to 0 (BBD exact), 1 "
148  "(inexact BBD with LU),"
149  << " 2 (inexact BBD with AMG) or 3 (exact BD).";
150  throw OomphLibError(
151  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
152  }
153 #endif
154 
155  // create the preconditioners
156  bool use_amg = true;
157  bool retain_all_blocks = false;
158  switch (Preconditioner_type)
159  {
160  // Exact BBD
161  case 0:
162 
163  retain_all_blocks = false;
165  new ExactSubBiharmonicPreconditioner(this, retain_all_blocks);
167 
168  oomph_info << "Using exact BBD\n";
169  break;
170 
171  // Inexact BBD with LU
172  case 1:
173 
174  use_amg = false;
176  new InexactSubBiharmonicPreconditioner(this, use_amg);
178  oomph_info << "Using inexact BBD with LU\n";
179  break;
180 
181 
182  // Inexact BBD with AMG
183  case 2:
184 
185  use_amg = true;
187  new InexactSubBiharmonicPreconditioner(this, use_amg);
189  oomph_info << "Using inexact BBD with AMG\n";
190  break;
191 
192  /// Exact BD
193  case 3:
194 
195  retain_all_blocks = true;
197  new ExactSubBiharmonicPreconditioner(this, retain_all_blocks);
199 
200  oomph_info << "Using exact BD\n";
201  break;
202 
203  default:
204 
205  throw OomphLibError("Wrong type of preconditioner.",
206  OOMPH_CURRENT_FUNCTION,
207  OOMPH_EXCEPTION_LOCATION);
208  }
209 
210 
211  // setup sub preconditioner pt 1
213 
214  // get the matrix ans setup sub preconditioner pt 2
215  CRDoubleMatrix* j_33_pt = new CRDoubleMatrix;
216  this->get_block(3, 3, *j_33_pt);
217  Sub_preconditioner_2_pt->setup(j_33_pt);
218  delete j_33_pt;
219  j_33_pt = 0;
220 
221  // if the block preconditioner has 5 block types setup the preconditioner
222  // for the 5th block diagonal block (Matrix is also diagonal hence a
223  // diagonal preconditioner is sufficient in the exact biharmonic
224  // preconditioner case as well)
225  if (this->nblock_types() == 5)
226  {
227  // get the matrix for block J_33
228  CRDoubleMatrix* j_44_pt = new CRDoubleMatrix;
229  this->get_block(4, 4, *j_44_pt);
230 
231  // setup the hijacked sub preconditioner
234  delete j_44_pt;
235  j_44_pt = 0;
236  }
237  }
238 
239 
240  //============================================================================
241  /// preconditioner solve for the biharmonic preconditioner
242  //============================================================================
244  DoubleVector& z)
245  {
246  // zero z
247  z.initialise(0.0);
248 
249  // solve sub preconditioner 1
251 
252  // solve sub preconditioner 2
253  DoubleVector block_r;
254  get_block_vector(3, r, block_r);
255  DoubleVector block_z;
257  return_block_vector(3, block_z, z);
258 
259  // solve the hijacked sub block preconditioner if required
260  if (this->nblock_types() == 5)
261  {
262  block_r.clear();
263  block_z.clear();
264  get_block_vector(4, r, block_r);
266  block_z);
267  return_block_vector(4, block_z, z);
268  }
269  }
270 
271  //============================================================================
272  /// setup for the exact sub biharmonic preconditioner
273  //============================================================================
275  {
276  // clean up memory first
277  this->clean_up_memory();
278 
279  // setup
280  this->block_setup();
281 
282  // Number of block types
283  unsigned n_block_types = this->nblock_types();
284 
285  // check for required number of blocks
286 #ifdef PARANOID
287  if (n_block_types != 3)
288  {
289  std::ostringstream error_message;
290  error_message
291  << "This preconditioner requires 3 block types.\n"
292  << "It is sub preconditioner for the BiharmonicPreconditioner.\n";
293  throw OomphLibError(
294  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
295  }
296 #endif
297 
298  // Data type indicating which blocks from the preconditioner matrix we want
299  VectorMatrix<BlockSelector> required_blocks(n_block_types, n_block_types);
300 
301  // boolean indicating if we want the block or not, stored for readability.
302  // Initially this is set to true for all blocks. Later we select which
303  // blocks we do not want.
304  const bool want_block = true;
305  for (unsigned b_i = 0; b_i < n_block_types; b_i++)
306  {
307  for (unsigned b_j = 0; b_j < n_block_types; b_j++)
308  {
309  required_blocks[b_i][b_j].select_block(b_i, b_j, want_block);
310  }
311  }
312 
313  // Which blocks do we not want?
314  if (!Retain_all_blocks)
315  {
316  required_blocks[1][2].do_not_want_block();
317  required_blocks[2][1].do_not_want_block();
318  }
319 
320  // Get the preconditioner matrix as defined by required_blocks
321  CRDoubleMatrix preconditioner_matrix =
322  this->get_concatenated_block(required_blocks);
323 
324  // setup the preconditioner
326  Sub_preconditioner_pt->setup(&preconditioner_matrix);
327 
328  // preconditioner_matrix will now go out of scope (and is destroyed).
329  }
330 
331 
332  //============================================================================
333  /// preconditioner solve for the exact sub biharmonic preconditioner
334  //============================================================================
336  const DoubleVector& r, DoubleVector& z)
337  {
338  // vectors for use within the sub preconditioner
339  DoubleVector sub_r;
340  DoubleVector sub_z;
341 
342  // get the sub r vector
344 
345  // solve the preconditioner
347 
348  // return the sub z vector to the master z vector
350  }
351 
352 
353  //============================================================================
354  /// setup for the inexact sub biharmonic preconditioner
355  //============================================================================
357  {
358  // clean up memory first
359  this->clean_up_memory();
360 
361  // setup
362  this->block_setup();
363 
364  // Number of block types
365  unsigned n_block_types = this->nblock_types();
366 
367  // paranoid check for number of blocks
368 #ifdef PARANOID
369  if (n_block_types != 3)
370  {
371  std::ostringstream error_message;
372  error_message
373  << "This preconditioner requires 3 block types.\n"
374  << "It is sub preconditioner for the BiharmonicPreconditioner.\n";
375  throw OomphLibError(
376  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
377  }
378 #endif
379 
380  // required blocks
381  DenseMatrix<bool> required_blocks(n_block_types, n_block_types, false);
382  required_blocks(0, 0) = true;
383  required_blocks(0, 1) = true;
384  required_blocks(1, 0) = true;
385  required_blocks(1, 1) = true;
386  required_blocks(0, 2) = true;
387  required_blocks(2, 0) = true;
388  required_blocks(2, 2) = true;
389 
390  // Matrix of block matrix pointers
391  Matrix_of_block_pointers.resize(n_block_types, n_block_types, 0);
392 
393  // get the blocks
394  this->get_blocks(required_blocks, Matrix_of_block_pointers);
395 
396  // lump the matrix J_11
400 
401  delete Matrix_of_block_pointers(1, 1);
402  Matrix_of_block_pointers(1, 1) = 0;
403 
404  // lump the matrix J_22
408  delete Matrix_of_block_pointers(2, 2);
409  Matrix_of_block_pointers(2, 2) = 0;
410 
411  // compute the schur complement
413 
414  // create the preconditioner for the S00 Schur complement linear system
415  if (Use_amg)
416  {
417 #ifdef OOMPH_HAS_HYPRE
418  // Use Hypre Boomer AMG
422 #else
423  std::ostringstream error_message;
424  error_message << "Request AMG solver but oomph-lib does not have HYPRE";
425  throw OomphLibError(
426  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
427 #endif
428  }
429  else
430  {
432  }
433 
434  // setup the preconditioner
436 
437  // clean up
438  delete S_00_pt;
439  S_00_pt = 0;
440  }
441 
442  //============================================================================
443  /// computes the schur complement for the inexact sub biharmonic
444  /// preconditioner
445  //============================================================================
447  {
448  // if required get pointers to the vector components of J01 and J10
449  int* J_01_row_start = 0;
450  int* J_01_column_index = 0;
451  double* J_01_value = 0;
452  int* J_10_row_start = 0;
453  int* J_10_column_index = 0;
454 
455  // J_01 matrix
456  J_01_row_start = Matrix_of_block_pointers(0, 1)->row_start();
457  J_01_column_index = Matrix_of_block_pointers(0, 1)->column_index();
458  J_01_value = Matrix_of_block_pointers(0, 1)->value();
459 
460  // J_10 matrix
461  J_10_row_start = Matrix_of_block_pointers(1, 0)->row_start();
462  J_10_column_index = Matrix_of_block_pointers(1, 0)->column_index();
463 
464  // if required get pointers to the vector components of J01 and J10
465  int* J_02_row_start = 0;
466  int* J_02_column_index = 0;
467  double* J_02_value = 0;
468  int* J_20_row_start = 0;
469  int* J_20_column_index = 0;
470 
471  // J_02 matrix
472  J_02_row_start = Matrix_of_block_pointers(0, 2)->row_start();
473  J_02_column_index = Matrix_of_block_pointers(0, 2)->column_index();
474  J_02_value = Matrix_of_block_pointers(0, 2)->value();
475 
476  // J_20 matrix
477  J_20_row_start = Matrix_of_block_pointers(2, 0)->row_start();
478  J_20_column_index = Matrix_of_block_pointers(2, 0)->column_index();
479 
480  // get the inverse lumped vector of J_11 if required
481  double* J_11_lumped_and_inverted = 0;
482  J_11_lumped_and_inverted =
483  Lumped_J_11_preconditioner_pt->inverse_lumped_vector_pt();
484 
485  // get the inverse lumped vector of J_22 if required
486  double* J_22_lumped_and_inverted = 0;
487  J_22_lumped_and_inverted =
488  Lumped_J_22_preconditioner_pt->inverse_lumped_vector_pt();
489 
490  // size of J00 matrix (and S00 matrix)
491  unsigned J_00_nrow = Matrix_of_block_pointers(0, 0)->nrow();
492 
493  // vectors for the schur complement
494  Vector<int> S_00_row_start(J_00_nrow + 1);
495  Vector<int> S_00_column_index;
496  Vector<double> S_00_value;
497 
498  // number of elements in the x-dimension of the mesh
499  unsigned n_element_x =
501  dynamic_cast<BiharmonicPreconditioner*>(
503  ->bulk_element_mesh_pt())
504  ->nelement_in_dim(0);
505 
506  // nnz in schur complement (initialised to zero)
507  unsigned S_00_nnz = 0;
508 
509  // loop over columns of schur complement matrix
510  for (unsigned i = 0; i < J_00_nrow; i++)
511  {
512  // set column_start
513  S_00_row_start[i] = S_00_nnz;
514 
515  // loop over rows in schur complement matrix
516  // the schur complement matrix has 5 well defined bands thus we only
517  // perform matrix-matrix multiplication for these bands
518  //
519  // where the diagonal band is 0 :
520  //
521  // band 1 : -2*n_element_x +/- 5
522  // 2 : -n_element_x +/- 3
523  // 3 : 0 +/- 3
524  // 4 : n_element_x +/- 3
525  // 5 : 2*n_element_x +/- 5
526  //
527  // regardless of the type or combination of boundary conditions applied
528 
529  // Vector for postion of the bands in S_00
530  Vector<std::pair<int, int>> band_position(5);
531 
532  // compute the minimum and maximum positions of each band in terms of
533  // row number for column j
534  // note : static_cast used because max and min don't work on unsigned
535  band_position[0].first =
536  std::max(0, static_cast<int>(i - n_element_x * 2 - 5));
537  band_position[0].second =
538  std::max(0,
539  std::min(static_cast<int>(J_00_nrow - 1),
540  static_cast<int>(i - n_element_x * 2 + 5)));
541  band_position[1].first =
542  std::max(band_position[0].second + 1,
543  std::max(0, static_cast<int>(i - n_element_x - 3)));
544  band_position[1].second =
545  std::max(0,
546  std::min(static_cast<int>(J_00_nrow - 1),
547  static_cast<int>(i - n_element_x + 3)));
548  band_position[2].first = std::max(band_position[1].second + 1,
549  std::max(0, static_cast<int>(i - 3)));
550  band_position[2].second = std::max(
551  0, std::min(static_cast<int>(J_00_nrow - 1), static_cast<int>(i + 3)));
552  band_position[3].first =
553  std::max(band_position[2].second + 1,
554  std::max(0, static_cast<int>(i + n_element_x - 3)));
555  band_position[3].second =
556  std::max(0,
557  std::min(static_cast<int>(J_00_nrow - 1),
558  static_cast<int>(i + n_element_x + 3)));
559  band_position[4].first =
560  std::max(band_position[3].second + 1,
561  std::max(0, static_cast<int>(i + n_element_x * 2 - 5)));
562  band_position[4].second =
563  std::max(0,
564  std::min(static_cast<int>(J_00_nrow - 1),
565  static_cast<int>(i + n_element_x * 2 + 5)));
566 
567  // number of bands
568  unsigned n_band = 5;
569 
570  // loop over the bands
571  for (unsigned b = 0; b < n_band; b++)
572  {
573  // loop over the rows in band b
574  for (unsigned j = static_cast<unsigned>(band_position[b].first);
575  j <= static_cast<unsigned>(band_position[b].second);
576  j++)
577  {
578  ;
579 
580  // temporary value for the computation of S00(i,j)
581  double temp_value = Matrix_of_block_pointers(0, 0)->operator()(i, j);
582 
583  // iterate through non-zero entries of column j of A_10
584  for (int k = J_01_row_start[i]; k < J_01_row_start[i + 1]; k++)
585  {
586  if (J_10_column_index[J_10_row_start[J_01_column_index[k]]] <=
587  static_cast<int>(j) &&
588  static_cast<int>(j) <=
589  J_10_column_index[J_10_row_start[J_01_column_index[k] + 1] -
590  1])
591  {
592  temp_value -= J_01_value[k] *
593  Matrix_of_block_pointers(1, 0)->operator()(
594  J_01_column_index[k], j) *
595  J_11_lumped_and_inverted[J_01_column_index[k]];
596  }
597  }
598 
599  // next compute contribution for A_02*lumped(A_22)'*A_20
600 
601  // iterate through non-zero entries of column j of A_10
602  for (int k = J_02_row_start[i]; k < J_02_row_start[i + 1]; k++)
603  {
604  if (J_20_column_index[J_20_row_start[J_02_column_index[k]]] <=
605  static_cast<int>(j) &&
606  static_cast<int>(j) <=
607  J_20_column_index[J_20_row_start[J_02_column_index[k] + 1] -
608  1])
609  {
610  temp_value -= J_02_value[k] *
611  Matrix_of_block_pointers(2, 0)->operator()(
612  J_02_column_index[k], j) *
613  J_22_lumped_and_inverted[J_02_column_index[k]];
614  }
615  }
616 
617  // add element to schur complement matrix S00
618  if (temp_value != 0.0)
619  {
620  S_00_nnz++;
621  S_00_value.push_back(temp_value);
622  S_00_column_index.push_back(j);
623  }
624  }
625  }
626  }
627 
628  // last entry of s00 column start
629  S_00_row_start[J_00_nrow] = S_00_nnz;
630 
631  // build the schur complement S00
633  J_00_nrow,
634  S_00_value,
635  S_00_column_index,
636  S_00_row_start);
637 
638  // replace block J01 with J01*lumped(J11)' (if J11 can be lumped)
639  unsigned J_01_nnz = Matrix_of_block_pointers(0, 1)->nnz();
640  for (unsigned i = 0; i < J_01_nnz; i++)
641  {
642  J_01_value[i] *= J_11_lumped_and_inverted[J_01_column_index[i]];
643  }
644 
645  // replace block J_02 with J_02*lumped(J_22)' (if J22 can be lumped)
646  unsigned J_02_nnz = Matrix_of_block_pointers(0, 2)->nnz();
647  for (unsigned i = 0; i < J_02_nnz; i++)
648  {
649  J_02_value[i] *= J_22_lumped_and_inverted[J_02_column_index[i]];
650  }
651  }
652 
653 
654  //============================================================================
655  /// preconditioner solve for the inexact sub biharmonic preconditioner
656  //============================================================================
658  const DoubleVector& r, DoubleVector& z)
659  {
660  // get the block vectors
661  Vector<DoubleVector> block_r(3);
662  get_block_vectors(r, block_r);
663 
664  // r_0 = r_0 - J_01 * lumped(J_11)'*r_1 - J_02 * lumped(J_22)'*r_2
665  // Remember that J_01 has already been premultiplied by lumped(J_11)
666  DoubleVector temp;
667  Matrix_of_block_pointers(0, 1)->multiply(block_r[1], temp);
668  block_r[0] -= temp;
669  temp.clear();
670  Matrix_of_block_pointers(0, 2)->multiply(block_r[2], temp);
671  block_r[0] -= temp;
672 
673  // apply the inexact preconditioner
674  temp.clear();
675  S_00_preconditioner_pt->preconditioner_solve(block_r[0], temp);
676  return_block_vector(0, temp, z);
677 
678  // solve: lumped(J_11) x_1 = r_1 - J_10 x_0 for x_1
679  // remember temp contains r_0 (...or z_0)
680  DoubleVector temp2;
681  Matrix_of_block_pointers(1, 0)->multiply(temp, temp2);
682  block_r[1] -= temp2;
683  DoubleVector z_1;
684  Lumped_J_11_preconditioner_pt->preconditioner_solve(block_r[1], z_1);
685  return_block_vector(1, z_1, z);
686 
687  // solve: lumped(J_22) x_2 = r_2 - J_20 x_0 for x_2
688  // remember temp contains r_0 (...or z_0)
689  temp2.clear();
690  Matrix_of_block_pointers(2, 0)->multiply(temp, temp2);
691  block_r[2] -= temp2;
692  DoubleVector z_2;
693  Lumped_J_22_preconditioner_pt->preconditioner_solve(block_r[2], z_2);
694  return_block_vector(2, z_2, z);
695  }
696 } // namespace oomph
cstr elem_len * i
Definition: cfortran.h:603
Biharmonic Preconditioner - for two dimensional problems.
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.
void setup()
Setup the preconditioner.
unsigned Preconditioner_type
preconditioner type
Mesh * Bulk_element_mesh_pt
the bulk element mesh pt
Preconditioner * Sub_preconditioner_1_pt
Exact Preconditioner Pointer.
void return_block_vector(const unsigned &n, const DoubleVector &b, DoubleVector &v) const
Takes the n-th block ordered vector, b, and copies its entries to the appropriate entries in the natu...
CRDoubleMatrix get_concatenated_block(const VectorMatrix< BlockSelector > &selected_block)
Returns a concatenation of the block matrices specified by the argument selected_block....
void get_blocks(DenseMatrix< bool > &required_blocks, DenseMatrix< CRDoubleMatrix * > &block_matrix_pt) const
Get all the block matrices required by the block preconditioner. Takes a pointer to a matrix of bools...
void get_block(const unsigned &i, const unsigned &j, CRDoubleMatrix &output_matrix, const bool &ignore_replacement_block=false) const
Put block (i,j) into output_matrix. This block accounts for any coarsening of dof types and any repla...
const LinearAlgebraDistribution * block_distribution_pt(const unsigned &b) const
Access function to the block distributions (const version).
void get_block_vector(const unsigned &n, const DoubleVector &v, DoubleVector &b) const
Takes the naturally ordered vector, v and returns the n-th block vector, b. Here n is the block numbe...
unsigned nblock_types() const
Return the number of block types.
void return_block_ordered_preconditioner_vector(const DoubleVector &w, DoubleVector &v) const
Takes the block ordered vector, w, and reorders it in natural order. Reordered vector is returned in ...
CRDoubleMatrix * matrix_pt() const
Access function to matrix_pt. If this is the master then cast the matrix pointer to MATRIX*,...
BlockPreconditioner< CRDoubleMatrix > * master_block_preconditioner_pt() const
Access function to the master block preconditioner pt.
void get_block_vectors(const Vector< unsigned > &block_vec_number, const DoubleVector &v, Vector< DoubleVector > &s) const
Takes the naturally ordered vector and rearranges it into a vector of sub vectors corresponding to th...
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...
void get_block_ordered_preconditioner_vector(const DoubleVector &v, DoubleVector &w)
Given the naturally ordered vector, v, return the vector rearranged in block order in w....
A class for compressed row matrices. This is a distributable object.
Definition: matrices.h:888
A vector in the mathematical sense, initially developed for linear algebra type applications....
Definition: double_vector.h:58
void initialise(const double &v)
initialise the whole vector with value v
void clear()
wipes the DoubleVector
Sub Biharmonic Preconditioner - an exact preconditioner for the 3x3 top left hand corner sub block ma...
void clean_up_memory()
delete the subsidiary preconditioner pointer
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.
A two dimensional Hermite bicubic element quadrilateral mesh for a topologically rectangular domain....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
Definition: hypre_solver.h:826
unsigned & amg_smoother_iterations()
Access function to AMG_smoother_iterations.
unsigned & amg_simple_smoother()
Access function to AMG_simple_smoother flag.
Definition: hypre_solver.h:983
double & amg_strength()
Access function to AMG_strength.
unsigned & amg_iterations()
Return function for Max_iter.
Definition: hypre_solver.h:970
unsigned & hypre_method()
Access function to Hypre_method flag – specified via enumeration.
Definition: hypre_solver.h:945
unsigned & amg_coarsening()
Access function to AMG_coarsening flag.
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...
MatrixBasedLumpedPreconditioner< CRDoubleMatrix > * Lumped_J_22_preconditioner_pt
Preconditioner for storing the lumped J_22 matrix.
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 diagonal preconditioner.
An OomphLibError object which should be thrown when an run-time error is encountered....
void setup(DoubleMatrixBase *matrix_pt)
Setup the preconditioner: store the matrix pointer and the communicator pointer then call preconditio...
virtual void preconditioner_solve(const DoubleVector &r, DoubleVector &z)=0
Apply the preconditioner. Pure virtual generic interface function. This method should apply the preco...
An interface to allow SuperLU to be used as an (exact) Preconditioner.
VectorMatrix is a generalised, STL-map-based, matrix based on a Vector of Vectors.
Definition: vector_matrix.h:79
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
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...