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