spacetime_navier_stokes_block_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-2024 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 // Header file
28 
29 /// /////////////////////////////////////////////////////////////////////////
30 /// /////////////////////////////////////////////////////////////////////////
31 /// /////////////////////////////////////////////////////////////////////////
32 
33 namespace oomph
34 {
35  /*
36  #ifdef OOMPH_HAS_HYPRE
37  //======start_of_HypreSubsidiaryPreconditionerHelper_namespace==============
38  /// Helper method for the block diagonal F block preconditioner to allow
39  /// hypre to be used as a subsidiary block preconditioner
40  //==========================================================================
41  namespace HypreSubsidiaryPreconditionerHelper
42  {
43  /// Assign the Hypre preconditioner pointer
44  Preconditioner* set_hypre_preconditioner()
45  {
46  return new HyprePreconditioner;
47  }
48  } // End of HypreSubsidiaryPreconditionerHelper
49  #endif
50  */
51 
52  //============================================================================
53  /// Setup for the block triangular preconditioner
54  //============================================================================
56  {
57  // For debugging...
58  bool doc_block_matrices = false;
59 
60  // Clean up any previously created data
62 
63  // If we're meant to build silently
64  if (this->Silent_preconditioner_setup == true)
65  {
66  // Store the output stream pointer
67  this->Stream_pt = oomph_info.stream_pt();
68 
69  // Now set the oomph_info stream pointer to the null stream to
70  // disable all possible output
72  }
73 
74 #ifdef PARANOID
75  // If the upcast was unsuccessful
76  if (dynamic_cast<CRDoubleMatrix*>(this->matrix_pt()) == 0)
77  {
78  // Allocate space for an error message
79  std::ostringstream error_message;
80 
81  // Create the error message
82  error_message << "NavierStokesSchurComplementPreconditioner only works "
83  << "with CRDoubleMatrix matrices" << std::endl;
84 
85  // Throw the error message
86  throw OomphLibError(
87  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
88  }
89 #endif
90 
91  // Start the timer for the block setup
92  double t_block_setup_start = TimingHelpers::timer();
93 
94  // Set up block look up schemes (done automatically in the
95  // BlockPreconditioner base class, based on the information
96  // provided in the block-preconditionable elements in the problem)
97 
98  // This preconditioner has two types of block:
99  // Type 0: Velocity - Corresponding to DOF type 0 & 1
100  // Type 1: Pressure - Corresponding to DOF type 2
101  unsigned n_dof_types = 0;
102 
103  // If this has been set up as a subsidiary preconditioner
105  {
106  // Get the number of dof types (will be told by the master preconditioner)
107  n_dof_types = this->ndof_types();
108 
109  // DRAIG: Delete
110  oomph_info << "Number of dof types: " << n_dof_types << std::endl;
111 
112  // Make sure there's only 3 dof types for now!
113  if (n_dof_types != 3)
114  {
115  // Allocate space for an error message
116  std::ostringstream error_message;
117 
118  // Create the error message
119  error_message << "Should only be 3 dof types! You have " << n_dof_types
120  << " dof types!" << std::endl;
121 
122  // Throw an error
123  throw OomphLibError(error_message.str(),
124  OOMPH_CURRENT_FUNCTION,
125  OOMPH_EXCEPTION_LOCATION);
126  }
127  }
128  else
129  {
130  // Throw an error
131  throw OomphLibError("Can only be used as a subsidiary preconditioner!",
132  OOMPH_CURRENT_FUNCTION,
133  OOMPH_EXCEPTION_LOCATION);
134  } // if (this->is_subsidiary_block_preconditioner())
135 
136  // Allocate storage for the dof-to-block mapping
137  Vector<unsigned> dof_to_block_map(n_dof_types);
138 
139  // The last dof type (the pressure) should have it's own block type
140  dof_to_block_map[n_dof_types - 1] = 1;
141 
142  // Call the generic block setup function
143  this->block_setup(dof_to_block_map);
144 
145  // Stop the timer
146  double t_block_setup_finish = TimingHelpers::timer();
147 
148  // Calculate the time difference
149  double block_setup_time = t_block_setup_finish - t_block_setup_start;
150 
151  // Document the time taken
152  oomph_info << "Time for block_setup(...) [sec]: " << block_setup_time
153  << std::endl;
154 
155  // Get the number of block types
156  unsigned n_block_type = this->nblock_types();
157 
158  // How many block types are there? Should only be 2...
159  oomph_info << "Number of block types: " << n_block_type << std::endl;
160 
161  // Allocate storage for the momentum block
162  CRDoubleMatrix* f_pt = new CRDoubleMatrix;
163 
164  // Allocate storage for the divergence block
165  CRDoubleMatrix* d_pt = new CRDoubleMatrix;
166 
167  // Allocate storage for the gradient block
168  CRDoubleMatrix* g_pt = new CRDoubleMatrix;
169 
170  // Allocate storage for the pressure Poisson matrix
171  CRDoubleMatrix* p_matrix_pt = new CRDoubleMatrix;
172 
173  // Get the momentum block F
174  this->get_block(0, 0, *f_pt);
175 
176  // Get the divergence block, D
177  this->get_block(1, 0, *d_pt);
178 
179  // Get gradient matrix, G=D^T
180  this->get_block(0, 1, *g_pt);
181 
182  // Calculate the pressure Poisson-like matrix
183  d_pt->multiply(*g_pt, *p_matrix_pt);
184 
185  // If the user wants the block matrices themselves
186  if (doc_block_matrices)
187  {
188  // Filename suffix (defaults should be an empty string)
189  std::string suffix = "";
190 
191  // Create space for all of the blocks
192  DenseMatrix<CRDoubleMatrix*> block_matrix_pt(
193  n_block_type, n_block_type, 0);
194 
195  // Loop over the block rows
196  for (unsigned i = 0; i < n_block_type; i++)
197  {
198  // Loop over the block columns
199  for (unsigned j = 0; j < n_block_type; j++)
200  {
201  // Allocate space for the matrix
202  block_matrix_pt(i, j) = new CRDoubleMatrix;
203 
204  // Get the block
205  this->get_block(i, j, *block_matrix_pt(i, j));
206  }
207  } // for (unsigned i=0;i<n_block_type;i++)
208 
209  // Space for the full matrix
210  CRDoubleMatrix* j_pt = new CRDoubleMatrix;
211 
212  // Concatenate the sub-blocks to get the full matrix
213  CRDoubleMatrixHelpers::concatenate(block_matrix_pt, *j_pt);
214 
215  // Output the matrix in MATLAB form
216  j_pt->sparse_indexed_output_with_offset("j_" + suffix + ".csv");
217 
218  // Tell the user
219  oomph_info << "Done output of j_" + suffix + ".csv" << std::endl;
220 
221  // Loop over the block rows
222  for (unsigned i = 0; i < n_block_type; i++)
223  {
224  // Loop over the block columns
225  for (unsigned j = 0; j < n_block_type; j++)
226  {
227  // Clear up this block
228  delete block_matrix_pt(i, j);
229 
230  // Make it a null pointer
231  block_matrix_pt(i, j) = 0;
232  }
233  } // for (unsigned i=0;i<n_block_type;i++)
234 
235  // Now clear up j_pt
236  delete j_pt;
237 
238  // Make it a null pointer
239  j_pt = 0;
240 
241  // Output the matrix in MATLAB form
242  f_pt->sparse_indexed_output_with_offset("f_matrix" + suffix + ".dat");
243 
244  // Tell the user
245  oomph_info << "Done output of f_matrix" + suffix + ".dat" << std::endl;
246 
247  // Output the matrix in MATLAB form
248  d_pt->sparse_indexed_output_with_offset("d_matrix" + suffix + ".dat");
249 
250  // Tell the user
251  oomph_info << "Done output of d_matrix" + suffix + ".dat" << std::endl;
252 
253  // Output the matrix in MATLAB form
254  g_pt->sparse_indexed_output_with_offset("g_matrix" + suffix + ".dat");
255 
256  // Tell the user
257  oomph_info << "Done output of g_matrix" + suffix + ".dat" << std::endl;
258 
259  // Output the matrix in MATLAB form
260  p_matrix_pt->sparse_indexed_output_with_offset("p_matrix" + suffix +
261  ".dat");
262 
263  // Tell the user
264  oomph_info << "Done output of p_matrix" + suffix + ".dat" << std::endl;
265 
266  // Exit
267  exit(0);
268  }
269 
270  // Kill divergence matrix because we don't need it any more
271  delete d_pt;
272 
273  // Make it a null pointer
274  d_pt = 0;
275 
276  // Make a new matrix-vector product
277  F_mat_vec_pt = new MatrixVectorProduct;
278 
279  // Make a new matrix-vector product
280  G_mat_vec_pt = new MatrixVectorProduct;
281 
282  // Set the matrix-vector product up
284 
285  // Set the matrix-vector product up
287 
288  // Do we need to doc. the memory statistics?
290  {
291  // Initialise the memory usage variable
292  Memory_usage_in_bytes = 0.0;
293 
294  // How many rows are there in the momentum block?
295  unsigned n_row_f = f_pt->nrow();
296 
297  // How many nonzeros are there in the momentum block?
298  unsigned n_nnz_f = f_pt->nnz();
299 
300  // Update the memory usage variable. NOTE: This calculation is done by
301  // taking into account the storage scheme for a compressed row matrix
302  Memory_usage_in_bytes += ((2 * ((n_row_f + 1) * sizeof(int))) +
303  (n_nnz_f * (sizeof(int) + sizeof(double))));
304 
305  // How many rows are there in the gradient block?
306  unsigned n_row_g = g_pt->nrow();
307 
308  // How many nonzeros are there in the gradient block?
309  unsigned n_nnz_g = g_pt->nnz();
310 
311  // Update the memory usage variable
312  Memory_usage_in_bytes += ((2 * ((n_row_g + 1) * sizeof(int))) +
313  (n_nnz_g * (sizeof(int) + sizeof(double))));
314  }
315 
316  // Kill the gradient matrix because we don't need it any more
317  delete g_pt;
318 
319  // Make it a null pointer
320  g_pt = 0;
321 
322  // If the P preconditioner has not been set
323  if (P_preconditioner_pt == 0)
324  {
325  // Just use SuperLU
326  P_preconditioner_pt = new SuperLUPreconditioner;
327 
328  // Indicate that we're using the default preconditioner
330 
331  // Do we need to doc. the memory statistics?
333  {
334  // Add on the memory needed to store and calculate the LU factors of P
336  dynamic_cast<SuperLUPreconditioner*>(P_preconditioner_pt)
337  ->get_total_memory_needed_for_superlu();
338  }
339  } // if (P_preconditioner_pt==0)
340 
341  // Compute the LU decomposition of P
342  P_preconditioner_pt->setup(p_matrix_pt);
343 
344  // Now delete the actual matrix P, we don't need it anymore
345  delete p_matrix_pt;
346 
347  // Make it a null pointer
348  p_matrix_pt = 0;
349 
350  // If the F preconditioner has not been setup
351  if (F_preconditioner_pt == 0)
352  {
353  // Just use SuperLU
354  F_preconditioner_pt = new SuperLUPreconditioner;
355 
356  // Indicate that we're using the default preconditioner
358 
359  // Do we need to doc. the memory statistics?
361  {
362  // Add on the memory needed to store and calculate the LU factors of P
364  dynamic_cast<SuperLUPreconditioner*>(F_preconditioner_pt)
365  ->get_total_memory_needed_for_superlu();
366  }
367  } // if (F_preconditioner_pt==0)
368 
369  // Compute the LU decomposition of F
370  F_preconditioner_pt->setup(f_pt);
371 
372  // Now delete the matrix, we don't need it anymore
373  delete f_pt;
374 
375  // Make it a null pointer
376  f_pt = 0;
377 
378  // Remember that the preconditioner has been setup so the stored
379  // information can be wiped when we come here next...
381 
382  // If we're meant to build silently, reassign the oomph stream pointer
383  if (this->Silent_preconditioner_setup == true)
384  {
385  // Store the output stream pointer
386  oomph_info.stream_pt() = this->Stream_pt;
387 
388  // Reset our own stream pointer
389  this->Stream_pt = 0;
390  }
391  } // End of setup
392 
393  //=============================================================================
394  /// Preconditioner solve for the block triangular preconditioner
395  //=============================================================================
397  const DoubleVector& r, DoubleVector& z)
398  {
399  // If z is not set up then give it the same distribution
400  if (!z.distribution_pt()->built())
401  {
402  // Build z with the same distribution as r
403  z.build(r.distribution_pt(), 0.0);
404  }
405 
406  // -----------------------------------------------------------------------
407  // Step 1 - apply approximate Schur inverse to pressure unknowns (block 1)
408  // -----------------------------------------------------------------------
409  // First working vectors
410  DoubleVector temp_vec;
411 
412  // Second working vector
413  DoubleVector another_temp_vec;
414 
415  // Copy pressure values from residual vector to temp_vec:
416  // Loop over all entries in the global vector (this one
417  // includes velocity and pressure dofs in some random fashion)
418  this->get_block_vector(1, r, temp_vec);
419 
420  // Compute the vector P^{-1}r_p
421  P_preconditioner_pt->preconditioner_solve(temp_vec, another_temp_vec);
422 
423  // Now clear the vector storing r_p
424  temp_vec.clear();
425 
426  // Now compute G(P^{-1}r_p)
427  G_mat_vec_pt->multiply(another_temp_vec, temp_vec);
428 
429  // Clear the vector storing P^{-1}r_p
430  another_temp_vec.clear();
431 
432  // Now update to get F(GP^{-1}r_p)
433  F_mat_vec_pt->multiply(temp_vec, another_temp_vec);
434 
435  // Clear the vector storing GP^{-1}r_p
436  temp_vec.clear();
437 
438  // Now compute D(FGP^{-1}r_p)=G^{T}(FGP^{-1}r_p)
439  G_mat_vec_pt->multiply_transpose(another_temp_vec, temp_vec);
440 
441  // Clear the vector storing FGP^{-1}r_p
442  another_temp_vec.clear();
443 
444  // Finally, compute the full approximation to -z_p;
445  // -z_p=P^{-1}(DFG)P^{-1}r_p NOTE: We'll move the sign over in the next step
446  P_preconditioner_pt->preconditioner_solve(temp_vec, another_temp_vec);
447 
448  // Rebuild temp_vec with another_temp_vec's distribution with all
449  // entries initialised to zero
450  temp_vec.build(another_temp_vec.distribution_pt(), 0.0);
451 
452  // Now fix the sign and store the result in temp_vec
453  temp_vec -= another_temp_vec;
454 
455  // Now copy temp_vec (i.e. z_p) back into the global vector z
456  return_block_vector(1, temp_vec, z);
457 
458  // ------------------------------------------------------------
459  // Step 2 - apply preconditioner to velocity unknowns (block 0)
460  // ------------------------------------------------------------
461  // Clear the vector storing z_p=-P^{-1}(DFG)P^{-1}r_p
462  temp_vec.clear();
463 
464  // Recall that another_temp_vec (computed above) contains the negative
465  // of the solution of the Schur complement system, -z_p. Multiply by G
466  // and store the result in temp_vec (vector resizes itself).
467  G_mat_vec_pt->multiply(another_temp_vec, temp_vec);
468 
469  // Now clear the vector storing P^{-1}(DFG)P^{-1}r_p
470  another_temp_vec.clear();
471 
472  // Get the residual vector entries associated with the velocities, r_u
473  get_block_vector(0, r, another_temp_vec);
474 
475  // Update the result to get r_u-Gz_p
476  another_temp_vec += temp_vec;
477 
478  // Compute the solution z_u which satisfies z_u=F^{-1}(r_u-Gz_p)
479  F_preconditioner_pt->preconditioner_solve(another_temp_vec, temp_vec);
480 
481  // Now store the result in the global solution vector
482  return_block_vector(0, temp_vec, z);
483  } // End of preconditioner_solve
484 
485 
486  //============================================================================
487  /// Setup for the GMRES block preconditioner
488  //============================================================================
490  {
491  // Clean up any previously created data
492  clean_up_memory();
493 
494  // If we're meant to build silently
495  if (this->Silent_preconditioner_setup == true)
496  {
497  // Store the output stream pointer
498  this->Stream_pt = oomph_info.stream_pt();
499 
500  // Now set the oomph_info stream pointer to the null stream to
501  // disable all possible output
503  }
504 
505  // Start the timer for the block setup
506  double t_block_setup_start = TimingHelpers::timer();
507 
508  // Set up block look up schemes (done automatically in the
509  // BlockPreconditioner base class, based on the information
510  // provided in the block-preconditionable elements in the problem)
511 
512  // This preconditioner has two types of block:
513  // Type 0: Velocity - Corresponding to DOF type 0 & 1
514  // Type 1: Pressure - Corresponding to DOF type 2
515  unsigned n_dof_types = 0;
516 
517  // If this has been set up as a subsidiary preconditioner
519  {
520  // Get the number of dof types (will be told by the master preconditioner)
521  n_dof_types = this->ndof_types();
522 
523  // Output the number of dof types
524  oomph_info << "Number of dof types: " << n_dof_types << std::endl;
525 
526  // Make sure there's only 3 dof types for now!
527  if (n_dof_types != 3)
528  {
529  // Allocate space for an error message
530  std::ostringstream error_message;
531 
532  // Create the error message
533  error_message << "Should only be 3 dof types! You have " << n_dof_types
534  << " dof types!" << std::endl;
535 
536  // Throw an error
537  throw OomphLibError(error_message.str(),
538  OOMPH_CURRENT_FUNCTION,
539  OOMPH_EXCEPTION_LOCATION);
540  }
541  }
542  // If it's being used as a master preconditioner
543  else
544  {
545  // Throw an error
546  throw OomphLibError("Currently only used as a subsidiary preconditioner!",
547  OOMPH_CURRENT_FUNCTION,
548  OOMPH_EXCEPTION_LOCATION);
549  } // if (this->is_subsidiary_block_preconditioner())
550 
551  // Aggregrate everything together since GMRES itself isn't going to take
552  // advantage of the block structure...
553  Vector<unsigned> dof_to_block_map(n_dof_types, 0);
554 
555  // Call the generic block setup function
556  this->block_setup(dof_to_block_map);
557 
558  // Stop the timer
559  double t_block_setup_finish = TimingHelpers::timer();
560 
561  // Calculate the time difference
562  double block_setup_time = t_block_setup_finish - t_block_setup_start;
563 
564  // Document the time taken
565  oomph_info << "Time for block_setup(...) [sec]: " << block_setup_time
566  << std::endl;
567 
568  // Get the number of block types
569  unsigned n_block_type = this->nblock_types();
570 
571  // How many block types are there? Should only be 2...
572  oomph_info << "Number of block types: " << n_block_type << std::endl;
573 
574  // Space for the concatenated block matrix
575  Matrix_pt = new CRDoubleMatrix;
576 
577  // Get the (one and only) block to be used by this block preconditioner
578  this->get_block(0, 0, *Matrix_pt);
579 
580  // Create a new instance of the NS subsidiary preconditioner
582  new SpaceTimeNavierStokesSubsidiaryPreconditioner;
583 
584  // Do we need to doc. the memory statistics?
586  {
587  // How many rows are there in this block?
588  unsigned n_row = Matrix_pt->nrow();
589 
590  // How many nonzeros are there in this block?
591  unsigned n_nnz = Matrix_pt->nnz();
592 
593  // Compute the memory usage. The only memory overhead here (for the
594  // setup phase) is storing the system matrix
595  Memory_usage_in_bytes = ((2 * ((n_row + 1) * sizeof(int))) +
596  (n_nnz * (sizeof(int) + sizeof(double))));
597 
598  // We might as well doc. the memory statistics of the Navier-Stokes
599  // subsidiary block preconditioner while we're at it...
601  }
602 
603  // Create a map to declare which of the 3 dof types in the GMRES
604  // subsidiary preconditioner correspond to the dof types in the NS
605  // subsidiary block preconditioner
606  Vector<unsigned> dof_map(n_dof_types);
607 
608  // Loop over the dof types associated with the GMRES subsidiary block
609  // preconditioner
610  for (unsigned i = 0; i < n_dof_types; i++)
611  {
612  // There is, essentially, an identity mapping between the GMRES subsidiary
613  // block preconditioner and the NSSP w.r.t. the dof types
614  dof_map[i] = i;
615  }
616 
617  // Now turn it into an "proper" subsidiary preconditioner
620 
621  // Setup: The NS subsidiary preconditioner is itself a block preconditioner
622  // so we need to pass it a pointer to full-size matrix!
624 
625  // Remember that the preconditioner has been setup so the stored
626  // information can be wiped when we come here next...
628 
629  // If we're meant to build silently, reassign the oomph stream pointer
630  if (this->Silent_preconditioner_setup == true)
631  {
632  // Store the output stream pointer
633  oomph_info.stream_pt() = this->Stream_pt;
634 
635  // Reset our own stream pointer
636  this->Stream_pt = 0;
637  }
638  } // End of setup
639 
640  //=============================================================================
641  /// Preconditioner solve for the GMRES block preconditioner
642  //=============================================================================
643  void GMRESBlockPreconditioner::preconditioner_solve(const DoubleVector& rhs,
644  DoubleVector& solution)
645  {
646  // Get the number of block types
647  unsigned n_block_types = this->nblock_types();
648 
649  // If there is more than one block type (there shouldn't be!)
650  if (n_block_types != 1)
651  {
652  // Throw an error
653  throw OomphLibError("Can only deal with one dof type!",
654  OOMPH_CURRENT_FUNCTION,
655  OOMPH_EXCEPTION_LOCATION);
656  }
657 
658  // Allocate space for the solution blocks
659  Vector<DoubleVector> block_z(n_block_types);
660 
661  // Allocate space for the residual blocks
662  Vector<DoubleVector> block_r(n_block_types);
663 
664  // Get the reordered RHS vector
665  this->get_block_vectors(rhs, block_r);
666 
667  // Get the reordered solution vector
668  this->get_block_vectors(solution, block_z);
669 
670  // Initialise the entries in the solution vector to zero
671  block_z[0].initialise(0.0);
672 
673  // Get number of dofs
674  unsigned n_dof = block_r[0].nrow();
675 
676  // Time solver
677  double t_start = TimingHelpers::timer();
678 
679  // Relative residual
680  double resid;
681 
682  // Iteration counter
683  unsigned iter = 1;
684 
685  // Initialise vectors
686  Vector<double> s(n_dof + 1, 0);
687  Vector<double> cs(n_dof + 1);
688  Vector<double> sn(n_dof + 1);
689  DoubleVector w(block_r[0].distribution_pt(), 0.0);
690 
691  // Storage for the time taken for all preconditioner applications
692  double t_prec_application_total = 0.0;
693 
694  // Storage for the RHS vector
695  DoubleVector r(block_r[0].distribution_pt(), 0.0);
696 
697  // If we're using LHS preconditioning then solve Mr=b-Jx for r (assumes x=0)
698  if (Preconditioner_LHS)
699  {
700  // Initialise timer
701  double t_prec_application_start = TimingHelpers::timer();
702 
703  // This is pretty inefficient but there is no other way to do this with
704  // block subsidiary preconditioners as far as I can tell because they
705  // demand to have the full r and z vectors...
706  DoubleVector block_z_with_size_of_full_z(rhs.distribution_pt(), 0.0);
707  DoubleVector r_updated(rhs.distribution_pt(), 0.0);
708 
709  // Reorder the entries in block_r[0] and put them back into the
710  // global-sized vector. The subsidiary preconditioner should only
711  // operate on the dofs that this preconditioner operates on so
712  // don't worry about all the other entries
713  this->return_block_vector(0, block_r[0], r_updated);
714 
715  // Solve on the global-sized vectors
717  r_updated, block_z_with_size_of_full_z);
718 
719  // Now grab the part of the solution associated with this preconditioner
720  this->get_block_vector(0, block_z_with_size_of_full_z, block_z[0]);
721 
722  // Doc time for setup of preconditioner
723  double t_prec_application_end = TimingHelpers::timer();
724 
725  // Update the preconditioner application time total
726  t_prec_application_total +=
727  (t_prec_application_end - t_prec_application_start);
728  }
729  // If we're using RHS preconditioning, do nothing to the actual RHS vector
730  else
731  {
732  // Copy the entries into r
733  r = block_r[0];
734  }
735 
736  // Calculate the initial residual norm (assumes x=0 initially)
737  double normb = r.norm();
738 
739  // Set beta (the initial residual)
740  double beta = normb;
741 
742  // Compute initial relative residual
743  if (normb == 0.0)
744  {
745  // To avoid dividing by zero, set normb to 1
746  normb = 1;
747  }
748 
749  // Now calculate the relative residual norm
750  resid = beta / normb;
751 
752  // If required, document convergence history to screen or file (if
753  // stream is open)
755  {
756  // If a file has not been provided to output to
757  if (!Output_file_stream.is_open())
758  {
759  // Output the initial relative residual norm to screen
760  oomph_info << 0 << " " << resid << std::endl;
761  }
762  // If there is a file provided to output to
763  else
764  {
765  // Output the initial relative residual norm to file
766  Output_file_stream << 0 << " " << resid << std::endl;
767  }
768  } // if (Doc_convergence_history)
769 
770  // If GMRES converges immediately
771  if (resid <= Tolerance)
772  {
773  // If the user wishes GMRES to output information about the solve to
774  // screen
775  if (Doc_time)
776  {
777  // Tell the user
778  oomph_info << "GMRES block preconditioner converged immediately. "
779  << "Normalised residual norm: " << resid << std::endl;
780  }
781 
782  // Now reorder the values in block_z[0] and put them back into the
783  // global solution vector
784  this->return_block_vector(0, block_z[0], solution);
785 
786  // Doc time for solver
787  double t_end = TimingHelpers::timer();
788  Solution_time = t_end - t_start;
789 
790  // If the user wishes GMRES to output information about the solve to
791  // screen
792  if (Doc_time)
793  {
794  // Tell the user
795  oomph_info << "Time for all preconditioner applications [sec]: "
796  << t_prec_application_total
797  << "\nTotal time for solve with GMRES block preconditioner "
798  << "[sec]: " << Solution_time << std::endl;
799  }
800 
801  // We're finished; end here!
802  return;
803  } // if (resid<=Tolerance)
804 
805  // Initialise vector of orthogonal basis vectors, v
806  Vector<DoubleVector> v(n_dof + 1);
807 
808  // Initialise the upper Hessenberg matrix H.
809  // NOTE: For implementation purposes the upper hessenberg matrix indices
810  // are swapped so the matrix is effectively transposed
811  Vector<Vector<double>> H(n_dof + 1);
812 
813  // Loop as many times as we're allowed
814  while (iter <= Max_iter)
815  {
816  // Copy the entries of the residual vector, r, into v[0]
817  v[0] = r;
818 
819  // Now update the zeroth basis vector, v[0], to have values r/beta
820  v[0] /= beta;
821 
822  // Storage for ...?
823  s[0] = beta;
824 
825  // Inner iteration counter for restarted version (we don't actually use
826  // a restart in this implementation of GMRES...)
827  unsigned iter_restart;
828 
829  // Perform iterations
830  for (iter_restart = 0; iter_restart < n_dof && iter <= Max_iter;
831  iter_restart++, iter++)
832  {
833  // Resize next column of the upper Hessenberg matrix
834  H[iter_restart].resize(iter_restart + 2);
835 
836  // Solve for w inside the scope of these braces
837  {
838  // Temporary vector
839  DoubleVector temp(block_r[0].distribution_pt(), 0.0);
840 
841  // If we're using LHS preconditioning solve Mw=Jv[i] for w
842  if (Preconditioner_LHS)
843  {
844  // Calculate temp=Jv[i]
845  Matrix_pt->multiply(v[iter_restart], temp);
846 
847  // Get the current time
848  double t_prec_application_start = TimingHelpers::timer();
849 
850  // This is pretty inefficient but there is no other way to do this
851  // with block sub preconditioners as far as I can tell because they
852  // demand to have the full r and z vectors...
853  DoubleVector block_z_with_size_of_full_z(rhs.distribution_pt(),
854  0.0);
855  DoubleVector r_updated(rhs.distribution_pt(), 0.0);
856 
857  // Put the values in temp into the global-sized vector
858  this->return_block_vector(0, temp, r_updated);
859 
860  // Apply the NSSP to the global-sized vectors (this will extract
861  // the appropriate sub-vectors and solve on the reduced system)
863  r_updated, block_z_with_size_of_full_z);
864 
865  // Now grab the part of the solution associated with this
866  // preconditioner
867  this->get_block_vector(0, block_z_with_size_of_full_z, w);
868 
869  // End timer
870  double t_prec_application_end = TimingHelpers::timer();
871 
872  // Update the preconditioner application time total
873  t_prec_application_total +=
874  (t_prec_application_end - t_prec_application_start);
875  }
876  // If we're using RHS preconditioning solve
877  else
878  {
879  // Initialise timer
880  double t_prec_application_start = TimingHelpers::timer();
881 
882  // This is pretty inefficient but there is no other way to do this
883  // with block sub preconditioners as far as I can tell because they
884  // demand to have the full r and z vectors...
885  DoubleVector block_z_with_size_of_full_z(rhs.distribution_pt());
886  DoubleVector r_updated(rhs.distribution_pt());
887 
888  // Put the values in v[iter_restart] back into the global-sized
889  // vector
890  this->return_block_vector(0, v[iter_restart], r_updated);
891 
892  // Solve on the global-sized vectors
894  r_updated, block_z_with_size_of_full_z);
895 
896  // Now grab the part of the solution associated with this
897  // preconditioner
898  this->get_block_vector(0, block_z_with_size_of_full_z, temp);
899 
900  // Doc time for setup of preconditioner
901  double t_prec_application_end = TimingHelpers::timer();
902 
903  // Update the preconditioner application time total
904  t_prec_application_total +=
905  (t_prec_application_end - t_prec_application_start);
906 
907  // Calculate w=Jv_m where v_m=M^{-1}v
908  Matrix_pt->multiply(temp, w);
909  }
910  } // End of solve for w
911 
912  // Get a pointer to the values in w
913  double* w_pt = w.values_pt();
914 
915  // Loop over the columns of the (transposed) Hessenberg matrix
916  for (unsigned k = 0; k <= iter_restart; k++)
917  {
918  // Initialise the entry in the upper Hessenberg matrix
919  H[iter_restart][k] = 0.0;
920 
921  // Get the pointer to the values of the k-th basis vector
922  double* vk_pt = v[k].values_pt();
923 
924  // Update H with the dot product of w and v[k]
925  H[iter_restart][k] += w.dot(v[k]);
926 
927  // Loop over the entries in w and v[k] again
928  for (unsigned i = 0; i < n_dof; i++)
929  {
930  // Now update the values in w
931  w_pt[i] -= H[iter_restart][k] * vk_pt[i];
932  }
933  } // for (unsigned k=0;k<=iter_restart;k++)
934 
935  // Calculate the (iter_restart+1,iter_restart)-th entry in the upper
936  // Hessenberg matrix (remember, H is actually the transpose of the
937  // proper upper Hessenberg matrix)
938  H[iter_restart][iter_restart + 1] = w.norm();
939 
940  // Build the (iter_restart+1)-th basis vector by copying w
941  v[iter_restart + 1] = w;
942 
943  // Now rescale the basis vector
944  v[iter_restart + 1] /= H[iter_restart][iter_restart + 1];
945 
946  // Loop over the columns of the transposed upper Hessenberg matrix
947  for (unsigned k = 0; k < iter_restart; k++)
948  {
949  // Apply a Givens rotation to entries of H
951  H[iter_restart][k], H[iter_restart][k + 1], cs[k], sn[k]);
952  }
953 
954  // Now generate the cs and sn entries for a plane rotation
955  generate_plane_rotation(H[iter_restart][iter_restart],
956  H[iter_restart][iter_restart + 1],
957  cs[iter_restart],
958  sn[iter_restart]);
959 
960  // Use the newly generated entries in cs and sn to apply a Givens
961  // rotation to those same entries
962  apply_plane_rotation(H[iter_restart][iter_restart],
963  H[iter_restart][iter_restart + 1],
964  cs[iter_restart],
965  sn[iter_restart]);
966 
967  // Now apply the Givens rotation to the entries in s
968  apply_plane_rotation(s[iter_restart],
969  s[iter_restart + 1],
970  cs[iter_restart],
971  sn[iter_restart]);
972 
973  // Compute the current residual
974  beta = std::fabs(s[iter_restart + 1]);
975 
976  // Compute the relative residual norm
977  resid = beta / normb;
978 
979  // If required, document the convergence history to screen or file
981  {
982  // If an output file has not been provided
983  if (!Output_file_stream.is_open())
984  {
985  // Output the convergence history to screen
986  oomph_info << iter << " " << resid << std::endl;
987  }
988  // An output file has been provided so output to that
989  else
990  {
991  // Output the convergence history to file
992  Output_file_stream << iter << " " << resid << std::endl;
993  }
994  } // if (Doc_convergence_history)
995 
996  // If the required tolerance was met
997  if (resid < Tolerance)
998  {
999  // Storage for the global-sized solution vector. Strictly speaking
1000  // we could actually just use the vector, solution but this can be
1001  // done after we know we've implemented everything correctly...
1002  DoubleVector block_z_with_size_of_full_z(rhs.distribution_pt(), 0.0);
1003 
1004  // Take the current solution, reorder the entries appropriately
1005  // and stick them into the global sized vector
1006  this->return_block_vector(0, block_z[0], block_z_with_size_of_full_z);
1007 
1008  // This will update block_z[0] and won't touch the global-size vector
1009  // block_x_with_size_of_full_x. We're in charge of reordering the
1010  // entries and putting it in solution
1011  update(
1012  iter_restart, H, s, v, block_z_with_size_of_full_z, block_z[0]);
1013 
1014  // Now go into block_z[0], reorder the entries in the appropriate
1015  // manner and put them into the vector, solution
1016  this->return_block_vector(0, block_z[0], solution);
1017 
1018  // If the user wishes GMRES to document the convergence
1019  if (Doc_time)
1020  {
1021  // Document the convergence to screen
1022  oomph_info
1023  << "\nGMRES block preconditioner converged (1). Normalised "
1024  << "residual norm: " << resid
1025  << "\nNumber of iterations to convergence: " << iter << "\n"
1026  << std::endl;
1027  }
1028 
1029  // Get the current time
1030  double t_end = TimingHelpers::timer();
1031 
1032  // Calculate the time difference, i.e. the time taken for the whole
1033  // solve
1034  Solution_time = t_end - t_start;
1035 
1036  // Storage the iteration count
1037  Iterations = iter;
1038 
1039  // If the user wishes GMRES to document the convergence
1040  if (Doc_time)
1041  {
1042  // Document the convergence to screen
1043  oomph_info
1044  << "Time for all preconditioner applications [sec]: "
1045  << t_prec_application_total
1046  << "\nTotal time for solve with GMRES block preconditioner "
1047  << "[sec]: " << Solution_time << std::endl;
1048  }
1049 
1050  // We're done now; finish here
1051  return;
1052  } // if (resid<Tolerance)
1053  } // for (iter_restart=0;iter_restart<n_dof&&iter<=Max_iter;...
1054 
1055  // Update
1056  if (iter_restart > 0)
1057  {
1058  // Storage for the global-sized solution vector
1059  DoubleVector block_z_with_size_of_full_z(rhs.distribution_pt(), 0.0);
1060 
1061  // Take the current solution, reorder the entries appropriately
1062  // and stick them into the global sized vector
1063  this->return_block_vector(0, block_z[0], block_z_with_size_of_full_z);
1064 
1065  // Update the (reordered block) solution vector
1066  update(
1067  (iter_restart - 1), H, s, v, block_z_with_size_of_full_z, block_z[0]);
1068 
1069  // Now go into block_z[0], reorder the entries in the appropriate manner
1070  // and put them into the vector, solution
1071  this->return_block_vector(0, block_z[0], solution);
1072  }
1073 
1074  // Solve Mr=(b-Jx) for r
1075  {
1076  // Temporary vector to store b-Jx
1077  DoubleVector temp(block_r[0].distribution_pt(), 0.0);
1078 
1079  // Calculate temp=b-Jx
1080  Matrix_pt->residual(block_z[0], block_r[0], temp);
1081 
1082  // If we're using LHS preconditioning
1083  if (Preconditioner_LHS)
1084  {
1085  // Initialise timer
1086  double t_prec_application_start = TimingHelpers::timer();
1087 
1088  // This is pretty inefficient but there is no other way to do this
1089  // with block sub preconditioners as far as I can tell because they
1090  // demand to have the full r and z vectors...
1091  DoubleVector block_z_with_size_of_full_z(rhs.distribution_pt());
1092  DoubleVector r_updated(rhs.distribution_pt());
1093 
1094  // Reorder the entries of temp and put them into the global-sized
1095  // vector
1096  this->return_block_vector(0, temp, r_updated);
1097 
1098  // Solve on the global-sized vectors
1100  r_updated, block_z_with_size_of_full_z);
1101 
1102  // Now grab the part of the solution associated with this
1103  // preconditioner
1104  this->get_block_vector(0, block_z_with_size_of_full_z, r);
1105 
1106  // Doc time for setup of preconditioner
1107  double t_prec_application_end = TimingHelpers::timer();
1108 
1109  // Update the preconditioner application time total
1110  t_prec_application_total +=
1111  (t_prec_application_end - t_prec_application_start);
1112  }
1113  } // End of solve Mr=(b-Jx) for r
1114 
1115  // Compute the current residual norm
1116  beta = r.norm();
1117 
1118  // if relative residual within tolerance
1119  resid = beta / normb;
1120  if (resid < Tolerance)
1121  {
1122  // Now store the result in the global solution vector
1123  this->return_block_vector(0, block_z[0], solution);
1124 
1125  if (Doc_time)
1126  {
1127  oomph_info
1128  << "\nGMRES block preconditioner converged (2). Normalised "
1129  << "residual norm: " << resid
1130  << "\nNumber of iterations to convergence: " << iter << "\n"
1131  << std::endl;
1132  }
1133 
1134  // Doc time for solver
1135  double t_end = TimingHelpers::timer();
1136  Solution_time = t_end - t_start;
1137 
1138  Iterations = iter;
1139 
1140  if (Doc_time)
1141  {
1142  oomph_info
1143  << "Time for all preconditioner applications [sec]: "
1144  << t_prec_application_total
1145  << "\nTotal time for solve with GMRES block preconditioner "
1146  << "[sec]: " << Solution_time << std::endl;
1147  }
1148  return;
1149  }
1150  }
1151 
1152 
1153  // otherwise GMRES failed convergence
1154  oomph_info << "\nGMRES block preconditioner did not converge to required "
1155  << "tolerance! \nReturning with normalised residual norm: "
1156  << resid << "\nafter " << Max_iter << " iterations.\n"
1157  << std::endl;
1158 
1160  {
1161  std::ostringstream error_message_stream;
1162  error_message_stream << "Solver failed to converge and you requested an "
1163  << "error on convergence failures.";
1164  throw OomphLibError(error_message_stream.str(),
1165  OOMPH_EXCEPTION_LOCATION,
1166  OOMPH_CURRENT_FUNCTION);
1167  }
1168  } // End of solve_helper
1169 } // End of namespace oomph
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
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...
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...
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.
CRDoubleMatrix * matrix_pt() const
Access function to matrix_pt. If this is the master then cast the matrix pointer to MATRIX*,...
bool is_subsidiary_block_preconditioner() const
Return true if this preconditioner is a subsidiary preconditioner.
unsigned ndof_types() const
Return the total number of DOF types.
void turn_into_subsidiary_block_preconditioner(BlockPreconditioner< MATRIX > *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 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...
void setup_matrix_vector_product(MatrixVectorProduct *matvec_prod_pt, CRDoubleMatrix *block_pt, const Vector< unsigned > &block_col_indices)
Setup a matrix vector product. matvec_prod_pt is a pointer to the MatrixVectorProduct,...
virtual void block_setup()
Determine the size of the matrix blocks and setup the lookup schemes relating the global degrees of f...
void multiply(const DoubleVector &x, DoubleVector &soln) const
Multiply the matrix by the vector x: soln=Ax.
Definition: matrices.cc:1782
unsigned long nnz() const
Return the number of nonzero entries (the local nnz)
Definition: matrices.h:1096
unsigned long nrow() const
Return the number of rows of the matrix.
Definition: matrices.h:1002
LinearAlgebraDistribution * distribution_pt() const
access to the LinearAlgebraDistribution
virtual void residual(const DoubleVector &x, const DoubleVector &b, DoubleVector &residual_)
Find the residual, i.e. r=b-Ax the residual.
Definition: matrices.h:326
double Memory_usage_in_bytes
Storage for the memory usage of the solver if the flag above is set to true (in bytes)
void apply_plane_rotation(double &dx, double &dy, double &cs, double &sn)
Helper function: Apply plane rotation. This is done using the update:
virtual void clean_up_memory()
Clean up the memory (empty for now...)
SpaceTimeNavierStokesSubsidiaryPreconditioner * Navier_stokes_subsidiary_preconditioner_pt
Pointer to the preconditioner for the block matrix.
bool Compute_memory_statistics
Flag to indicate whether or not to record the memory statistics this preconditioner.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r.
bool Preconditioner_has_been_setup
Control flag is true if the preconditioner has been setup (used so we can wipe the data when the prec...
bool Preconditioner_LHS
boolean indicating use of left hand preconditioning (if true) or right hand preconditioning (if false...
void update(const unsigned &k, const Vector< Vector< double >> &H, const Vector< double > &s, const Vector< DoubleVector > &v, const DoubleVector &block_x_with_size_of_full_x, DoubleVector &x)
Helper function to update the result vector using the result, x=x_0+V_m*y.
void generate_plane_rotation(double &dx, double &dy, double &cs, double &sn)
Helper function: Generate a plane rotation. This is done by finding the values of (i....
double Tolerance
Convergence tolerance.
bool Throw_error_after_max_iter
Should we throw an error instead of just returning when we hit the max iterations?
unsigned Max_iter
Maximum number of iterations.
double Solution_time
linear solver solution time
std::ofstream Output_file_stream
Output file stream for convergence history.
bool Doc_convergence_history
Flag indicating if the convergence history is to be documented.
bool Doc_time
Boolean flag that indicates whether the time taken.
Definition: linear_solver.h:77
void multiply_transpose(const DoubleVector &x, DoubleVector &y) const
Apply the transpose of the operator to the vector x and return the result in the vector y.
void multiply(const DoubleVector &x, DoubleVector &y) const
Apply the operator to the vector x and return the result in the vector y.
std::ostream *& stream_pt()
Access function for the stream pointer.
void setup(DoubleMatrixBase *matrix_pt)
Setup the preconditioner: store the matrix pointer and the communicator pointer then call preconditio...
std::ostream * Stream_pt
Pointer to the output stream – defaults to std::cout.
virtual void preconditioner_solve(const DoubleVector &r, DoubleVector &z)=0
Apply the preconditioner. Pure virtual generic interface function. This method should apply the preco...
bool Silent_preconditioner_setup
Boolean to indicate whether or not the build should be done silently.
Preconditioner * P_preconditioner_pt
Pointer to the 'preconditioner' for the pressure matrix.
double Memory_usage_in_bytes
Storage for the memory usage of the solver if the flag above is set to true (in bytes)
bool Preconditioner_has_been_setup
Control flag is true if the preconditioner has been setup (used so we can wipe the data when the prec...
bool Compute_memory_statistics
Flag to indicate whether or not to record the memory statistics this preconditioner.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r.
bool Using_default_f_preconditioner
Flag indicating whether the default F preconditioner is used.
bool Using_default_p_preconditioner
Flag indicating whether the default P preconditioner is used.
Preconditioner * F_preconditioner_pt
Pointer to the 'preconditioner' for the F matrix.
void concatenate(const DenseMatrix< CRDoubleMatrix * > &matrix_pt, CRDoubleMatrix &result_matrix)
Concatenate CRDoubleMatrix matrices. The in matrices are concatenated such that the block structure o...
Definition: matrices.cc:4349
std::string string(const unsigned &i)
Return the i-th string or "" if the relevant string hasn't been defined.
double timer()
returns the time in seconds after some point in past
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
Nullstream oomph_nullstream
///////////////////////////////////////////////////////////////////////// ///////////////////////////...
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...