pseudo_elastic_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-2025 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
28
29namespace oomph
30{
31 //=============================================================================
32 /// Functions to create instances of optimal subsidiary operators for
33 /// the PseudoElasticPreconditioner
34 //=============================================================================
35 namespace Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper
36 {
37#ifdef OOMPH_HAS_HYPRE
38
39 /// AMG w/ GS smoothing for the augmented elastic subsidiary linear
40 /// systems
42 {
44 new HyprePreconditioner("Hypre for diagonal blocks in pseudo-solid");
45 hypre_preconditioner_pt->set_amg_iterations(2);
46 hypre_preconditioner_pt->amg_using_simple_smoothing();
48 {
49 // Jacobi in parallel
50 hypre_preconditioner_pt->amg_simple_smoother() = 0;
51 }
52 else
53 {
54 // Gauss Seidel in serial (was 3 actually...)
55 hypre_preconditioner_pt->amg_simple_smoother() = 1;
56 }
58 hypre_preconditioner_pt->amg_damping() = 1.0;
59 hypre_preconditioner_pt->amg_coarsening() = 6;
61 }
62
63
64 /// AMG w/ GS smoothing for the augmented elastic subsidiary linear
65 /// systems -- calls Hypre version to stay consistent with previous default
70
71#endif
72
73#ifdef OOMPH_HAS_TRILINOS
74
75 /// TrilinosML smoothing for the augmented elastic
76 /// subsidiary linear systems
82
83 /// CG with diagonal preconditioner for the lagrange multiplier
84 /// subsidiary linear systems.
86 {
91 // Note: This makes CG a proper "inner iteration" for
92 // which GMRES (may) no longer converge. We should really
93 // use FGMRES or GMRESR for this. However, here the solver
94 // is so good that it'll converge very quickly anyway
95 // so there isn't much to be gained by limiting the number
96 // of iterations...
97 // prec_pt->max_iter() = 4;
98 prec_pt->solver_pt()->solver_type() = TrilinosAztecOOSolver::CG;
99 prec_pt->solver_pt()->disable_doc_time();
100 return prec_pt;
101 }
102#endif
103 } // namespace Pseudo_Elastic_Preconditioner_Subsidiary_Operator_Helper
104
105
106 ///////////////////////////////////////////////////////////////////////////////
107 ///////////////////////////////////////////////////////////////////////////////
108 ///////////////////////////////////////////////////////////////////////////////
109
110
111 //=============================================================================
112 // Setup method for the PseudoElasticPreconditioner.
113 //=============================================================================
115 {
116 // clean
117 this->clean_up_memory();
118
119#ifdef PARANOID
120 // paranoid check that meshes have been set
121 if (Elastic_mesh_pt == 0)
122 {
123 std::ostringstream error_message;
124 error_message << "The elastic mesh must be set.";
125 throw OomphLibError(
127 }
129 {
130 std::ostringstream error_message;
131 error_message << "The Lagrange multiplier mesh must be set.";
132 throw OomphLibError(
134 }
135#endif
136
137 // set the meshes
138 this->set_nmesh(2);
139 this->set_mesh(0, Elastic_mesh_pt);
141
142 // Figure out number of dof types
143 unsigned n_solid_dof_types = 0;
144 unsigned n_dof_types = 0;
145
147 {
148 // get the number of solid dof types from the first element
149 n_solid_dof_types = this->ndof_types_in_mesh(0);
150
151 // get the total number of dof types
152 n_dof_types = n_solid_dof_types + this->ndof_types_in_mesh(1);
153 }
154 else
155 {
156 n_dof_types = this->ndof_types();
157 n_solid_dof_types = n_dof_types - (n_dof_types / 3);
158 }
159#ifdef PARANOID
160 if (n_dof_types % 3 != 0)
161 {
162 std::ostringstream error_message;
163 error_message << "This preconditioner requires DIM*3 types of DOF";
164 throw OomphLibError(
166 }
167#endif
168
169 // determine the dimension
170 Dim = n_dof_types / 3;
171
172#ifdef PARANOID
173 // Recast Jacobian matrix to CRDoubleMatrix
175
176 if (cr_matrix_pt == 0)
177 {
178 std::ostringstream error_message;
179 error_message << "FSIPreconditioner only works with"
180 << " CRDoubleMatrix matrices" << std::endl;
181 throw OomphLibError(
183 }
184#endif
185
186 // Setup the dof_list scheme for block_setup.
187 // The natural DOF types order is (in 3D):
188 // 0 1 2 3 4 5 6 7 8 <- natural DOF type index
189 // x y z xc yc zc lx ly lz <- DOf type
190 //
191 // We need to group the directional displacements together:
192 // 0 1 2 3 4 5 6 7 8 <- block index
193 // x xc y yc xz zc lx ly lz <- DOF type
194 //
195 // The mapping required is:
196 // 0 2 4 1 3 5 6 7 8
197 //
198 // In general we want:
199 //
200 // dof_to_block_map[DOF type] = block type
202 for (unsigned dim_i = 0; dim_i < Dim; dim_i++)
203 {
204 // bulk solid dof types
206
207 // constrained solid dof types
209
210 // lagr dof types
212 }
213
214 // Setup the block ordering. We have the following:
215 // Block types [0, 2*Dim) are the solid blocks.
216 // Block types [2*Dim, 3*Dim) are the Lagrange multiplier blocks.
217 //
218 // For d = 0, 1, 2,
219 // Bulk solid doftypes: 2*d
220 // Constrained solid doftypes: 2*d+1
221 // Lgr doftyoes: 2*Dim + d
222 this->block_setup(dof_list_for_block_setup);
223
224 // Create dof_list for subsidiary preconditioner. This will be used later
225 // in turn_into_subsidiary_block_preconditioner(...)
227
228 // The dof types are currently in the order (in 3D):
229 // 0 1 2 3 4 5 6 7 8 <- DOF index
230 // x y z xc yc zc lx ly lz <- DOF type
231 //
232 // We need to group the directional displacements together:
233 // x xc y yc xz zc
234 //
235 // The mapping required is:
236 // 0 3 1 4 2 5
237 //
238 // In general we want:
239 // dof_list[subsidiary DOF type] = master DOF type
240 for (unsigned dim_i = 0; dim_i < Dim; dim_i++)
241 {
242 // bulk solid dof
244
245 // constrained solid dof
247 }
248
249 // Get the solid blocks
252
253 for (unsigned row_i = 0; row_i < n_solid_dof_types; row_i++)
254 {
255 for (unsigned col_i = 0; col_i < n_solid_dof_types; col_i++)
256 {
259 }
260 }
261
262 // compute the scaling (if required)
264 {
266 }
267 else
268 {
269 Scaling = 1.0;
270 }
271
272 // Add the scaled identify matrix to the constrained solid blocks.
273 for (unsigned d = 0; d < Dim; d++)
274 {
275 // Where is the constrained block located?
276 unsigned block_i = 2 * d + 1;
277
278 // Data from the constrained block.
279 double* s_values = solid_matrix_pt(block_i, block_i)->value();
280 int* s_column_index = solid_matrix_pt(block_i, block_i)->column_index();
281 int* s_row_start = solid_matrix_pt(block_i, block_i)->row_start();
282 int s_nrow_local = solid_matrix_pt(block_i, block_i)->nrow_local();
283 int s_first_row = solid_matrix_pt(block_i, block_i)->first_row();
284
285 // Loop through the diagonal entries and add the scaling.
286 for (int i = 0; i < s_nrow_local; i++)
287 {
288 bool found = false;
289 for (int j = s_row_start[i]; j < s_row_start[i + 1] && !found; j++)
290 {
291 if (s_column_index[j] == i + s_first_row)
292 {
293 s_values[j] += Scaling;
294 found = true;
295 } // if
296 } // for row start
297
298 // Check if the diagonal entry is found.
299 if (!found)
300 {
301 std::ostringstream error_message;
302 error_message << "The diagonal entry for the constained block("
303 << block_i << "," << block_i << ")\n"
304 << "on local row " << i << " does not exist."
305 << std::endl;
306 throw OomphLibError(error_message.str(),
309 }
310 } // for nrow_local
311 } // for Dim
312
313
314 // setup the solid subsidiary preconditioner
315 // //////////////////////////////// this preconditioner uses the full S
316 // matrix
318 {
321
322 // Add mesh (not actually used since this only acts as a subsidiary
323 // preconditioner...
324 // s_prec_pt->add_mesh(Elastic_mesh_pt);
325
327 Vector<unsigned>(1, 0));
328
329 for (unsigned i = 0; i < n_solid_dof_types; i++)
330 {
332 }
333
334 s_prec_pt->turn_into_subsidiary_block_preconditioner(
336
338 {
339 s_prec_pt->set_subsidiary_preconditioner_function(
341 }
342
343 // Set the replacement blocks.
344 for (unsigned d = 0; d < Dim; d++)
345 {
346 // The dof-block located in the block-ordering.
347 unsigned block_i = 2 * d + 1;
348
349 // The dof-block located in the dof-ordering.
350 unsigned dof_block_i = Dim + d;
353 }
354
355 s_prec_pt->Preconditioner::setup(matrix_pt());
357 }
358 // otherwise it is a block based preconditioner
359 else
360 {
362
363 // set the block preconditioning method
364 switch (E_preconditioner_type)
365 {
367 {
369 }
370 break;
372 {
376 block_triangular_prec_pt->upper_triangular();
377
379 }
380 break;
382 {
386 block_triangular_prec_pt->lower_triangular();
387
389 }
390 break;
391 default:
392 {
393 std::ostringstream error_msg;
395 << "There is no such block based preconditioner.\n"
396 << "Candidates are:\n"
397 << "PseudoElasticPreconditioner::Block_diagonal_preconditioner\n"
398 << "PseudoElasticPreconditioner::Block_upper_triangular_"
399 "preconditioner\n"
400 << "PseudoElasticPreconditioner::Block_lower_triangular_"
401 "preconditioner\n";
402 throw OomphLibError(
404 }
405 break;
406 }
407
408
409 // Add mesh (not actually used since this only acts as a subsidiary
410 // preconditioner...
411 // s_prec_pt->add_mesh(Elastic_mesh_pt);
412
413 // The block to block map
415 Vector<unsigned>(2, 0));
417
418 unsigned tmp_index = 0;
419 for (unsigned d = 0; d < Dim; d++)
420 {
424 }
425
426 s_prec_pt->turn_into_subsidiary_block_preconditioner(
428
430 {
431 s_prec_pt->set_subsidiary_preconditioner_function(
433 }
434
435 // Set the replacement blocks.
436 for (unsigned d = 0; d < Dim; d++)
437 {
438 // The dof-block located in the block-ordering.
439 unsigned block_i = 2 * d + 1;
440
441 // Then dof-block located in the dof-ordering.
442 unsigned dof_block_i = Dim + d;
445 }
446
447 s_prec_pt->set_dof_to_block_map(s_prec_dof_to_block_map);
448 s_prec_pt->Preconditioner::setup(matrix_pt());
449
451 }
452
453 // No longer require the solid blocks
454 for (unsigned row_i = 0; row_i < n_solid_dof_types; row_i++)
455 {
456 for (unsigned col_i = 0; col_i < n_solid_dof_types; col_i++)
457 {
458 delete solid_matrix_pt(row_i, col_i);
460 }
461 }
462
463 // next setup the lagrange multiplier preconditioners
465 for (unsigned d = 0; d < Dim; d++)
466 {
468 this->get_block(2 * Dim + d, 2 * d + 1, *b_pt);
469
470 // if a non default preconditioner is specified create
471 // the preconditioners
473 {
475 (*Lagrange_multiplier_subsidiary_preconditioner_function_pt)();
476 }
477 // else use default superlu preconditioner
478 else
479 {
482 }
483
484 // and setup
486 delete b_pt;
487 b_pt = 0;
488 }
489 }
490
491 //=============================================================================
492 /// Apply the elastic subsidiary preconditioner.
493 //=============================================================================
495 const DoubleVector& r, DoubleVector& z)
496 {
497 // apply the solid preconditioner
499 }
500
501 //=============================================================================
502 /// Apply the lagrange multiplier subsidiary preconditioner.
503 //=============================================================================
505 const DoubleVector& r, DoubleVector& z)
506 {
507 // apply the lagrange multiplier preconditioner
508 for (unsigned d = 0; d < Dim; d++)
509 {
510 DoubleVector x;
511 this->get_block_vector(Dim * 2 + d, r, x);
512 DoubleVector y;
513 Lagrange_multiplier_preconditioner_pt[d]->preconditioner_solve(x, y);
514 Lagrange_multiplier_preconditioner_pt[d]->preconditioner_solve(y, x);
515 unsigned nrow_local = x.nrow_local();
516 double* x_pt = x.values_pt();
517 for (unsigned i = 0; i < nrow_local; i++)
518 {
519 x_pt[i] = x_pt[i] * Scaling;
520 }
521 this->return_block_vector(Dim * 2 + d, x, z);
522 }
523 }
524
525 //=============================================================================
526 /// Clears the memory.
527 //=============================================================================
529 {
530 // clean the block preconditioner base class memory
532
533 // delete the solid preconditioner
536
537 // delete the lagrange multiplier preconditioner pt
539 for (unsigned i = 0; i < sz; i++)
540 {
543 }
544 }
545
546 ///////////////////////////////////////////////////////////////////////////////
547 ///////////////////////////////////////////////////////////////////////////////
548 ///////////////////////////////////////////////////////////////////////////////
549
550 //=============================================================================
551 // Setup method for the PseudoElasticPreconditioner.
552 //=============================================================================
554 {
555 // clean
556 this->clean_up_memory();
557
558#ifdef PARANOID
559 // paranoid check that meshes have been set
560 if (Elastic_mesh_pt == 0)
561 {
562 std::ostringstream error_message;
563 error_message << "The elastic mesh must be set.";
564 throw OomphLibError(
566 }
568 {
569 std::ostringstream error_message;
570 error_message << "The Lagrange multiplier mesh must be set.";
571 throw OomphLibError(
573 }
574#endif
575
576 // set the mesh
577 unsigned n_solid_dof_types = 0;
578 unsigned n_dof_types = 0;
579 this->set_mesh(0, Elastic_mesh_pt);
582 {
583 // get the number of solid dof types from the first element
584 n_solid_dof_types = this->ndof_types_in_mesh(0);
585
586 // get the total number of dof types
587 n_dof_types = n_solid_dof_types + this->ndof_types_in_mesh(1);
588 }
589 else
590 {
591 n_dof_types = this->ndof_types();
592 n_solid_dof_types = n_dof_types - (n_dof_types / 3);
593 }
594#ifdef PARANOID
595 if (n_dof_types % 3 != 0)
596 {
597 std::ostringstream error_message;
598 error_message << "This preconditioner requires DIM*3 types of DOF";
599 throw OomphLibError(
601 }
602#endif
603
604 // determine the dimension
605 Dim = n_dof_types / 3;
606
607 // Recast Jacobian matrix to CRDoubleMatrix
609
610#ifdef PARANOID
611 if (cr_matrix_pt == 0)
612 {
613 std::ostringstream error_message;
614 error_message << "FSIPreconditioner only works with"
615 << " CRDoubleMatrix matrices" << std::endl;
616 throw OomphLibError(
618 }
619#endif
620
621 // Setup the block ordering.
622 this->block_setup();
623
624 // Create dof_list for subsidiary preconditioner. This will be used later
625 // in turn_into_subsidiary_block_preconditioner(...)
627 for (unsigned i = 0; i < n_solid_dof_types; i++)
628 {
630 }
631
632 // compute the scaling (if required)
634 {
636 for (unsigned i = 0; i < n_solid_dof_types; i++)
637 {
638 dof_list[i] = i;
639 }
640
644 Scaling = helper_pt->s_inf_norm();
645 delete helper_pt;
646 helper_pt = 0;
647 }
648 else
649 {
650 Scaling = 1.0;
651 }
652
653
654 // setup the solid subsidiary preconditioner
655 // //////////////////////////////// this preconditioner uses the full S
656 // matrix
658 {
662 for (unsigned i = 0; i < n_solid_dof_types; i++)
663 {
664 dof_list[i] = i;
665 }
666 s_prec_pt->turn_into_subsidiary_block_preconditioner(this, dof_list);
668 {
669 s_prec_pt->set_subsidiary_preconditioner_function(
671 }
672 s_prec_pt->scaling() = Scaling;
673 s_prec_pt->Preconditioner::setup(matrix_pt());
675 }
676
677 // otherwise it is a block based preconditioner
678 else
679 {
680 // create the preconditioner
684 for (unsigned i = 0; i < n_solid_dof_types; i++)
685 {
686 dof_list[i] = i;
687 }
688 s_prec_pt->turn_into_subsidiary_block_preconditioner(this, dof_list);
689
690 // set the subsidiary solve method
692 {
693 s_prec_pt->set_subsidiary_preconditioner_function(
695 }
696
697 // set the scaling
698 s_prec_pt->scaling() = Scaling;
699
700 // set the block preconditioning method
701 switch (E_preconditioner_type)
702 {
704 s_prec_pt->use_block_diagonal_approximation();
705 break;
707 s_prec_pt->use_upper_triangular_approximation();
708 break;
710 s_prec_pt->use_lower_triangular_approximation();
711 break;
712 default:
713 break;
714 }
715
716 // setup
717 s_prec_pt->Preconditioner::setup(matrix_pt());
719 }
720
721 // next setup the lagrange multiplier preconditioners
723 for (unsigned d = 0; d < Dim; d++)
724 {
726 this->get_block(2 * Dim + d, Dim + d, *b_pt);
727
728 // if a non default preconditioner is specified create
729 // the preconditioners
731 {
733 (*Lagrange_multiplier_subsidiary_preconditioner_function_pt)();
734 }
735
736 // else use default superlu preconditioner
737 else
738 {
741 }
742
743 // and setup
745 delete b_pt;
746 b_pt = 0;
747 }
748 }
749
750 //=============================================================================
751 /// Apply the elastic subsidiary preconditioner.
752 //=============================================================================
754 const DoubleVector& r, DoubleVector& z)
755 {
756 // apply the solid preconditioner
758 }
759
760
761 //=============================================================================
762 /// Apply the lagrange multiplier subsidiary preconditioner.
763 //=============================================================================
765 const DoubleVector& r, DoubleVector& z)
766 {
767 // apply the lagrange multiplier preconditioner
768 for (unsigned d = 0; d < Dim; d++)
769 {
770 DoubleVector x;
771 this->get_block_vector(Dim * 2 + d, r, x);
772 DoubleVector y;
773 Lagrange_multiplier_preconditioner_pt[d]->preconditioner_solve(x, y);
774 Lagrange_multiplier_preconditioner_pt[d]->preconditioner_solve(y, x);
775 unsigned nrow_local = x.nrow_local();
776 double* x_pt = x.values_pt();
777 for (unsigned i = 0; i < nrow_local; i++)
778 {
779 x_pt[i] = x_pt[i] * Scaling;
780 }
781 this->return_block_vector(Dim * 2 + d, x, z);
782 }
783 }
784
785 //=============================================================================
786 /// Clears the memory.
787 //=============================================================================
789 {
790 // clean the block preconditioner base class memory
792
793 // delete the solid preconditioner
796
797 // delete the lagrange multiplier preconditioner pt
799 for (unsigned i = 0; i < sz; i++)
800 {
803 }
804 }
805
806
807 ///////////////////////////////////////////////////////////////////////////////
808 ///////////////////////////////////////////////////////////////////////////////
809 ///////////////////////////////////////////////////////////////////////////////
810
811
812 //=============================================================================
813 /// Setup the preconditioner
814 //=============================================================================
816 {
817 // clean memory
818 this->clean_up_memory();
819
820#ifdef PARANOID
821 // paranoid check that this preconditioner has an even number of DOF types
822 if (this->ndof_types() % 2 != 0)
823 {
824 std::ostringstream error_message;
825 error_message
826 << "This SUBSIDIARY preconditioner requires an even number of "
827 << "types of DOF";
828 throw OomphLibError(
830 }
831#endif
832
833 // assemble dof_to_block_map
834 unsigned ndof_types = this->ndof_types();
836 for (unsigned i = ndof_types / 2; i < ndof_types; i++)
837 {
838 dof_to_block_map[i] = 1;
839 }
840
841 this->block_setup(dof_to_block_map);
842
843 // get block 11
845 this->get_block(1, 1, *s11_pt);
846
847 // add the scaled identity matrix to block 11
848 double* s11_values = s11_pt->value();
849 int* s11_column_index = s11_pt->column_index();
850 int* s11_row_start = s11_pt->row_start();
851 int s11_nrow_local = s11_pt->nrow_local();
852 int s11_first_row = s11_pt->first_row();
853 for (int i = 0; i < s11_nrow_local; i++)
854 {
855 bool found = false;
856 for (int j = s11_row_start[i]; j < s11_row_start[i + 1] && !found; j++)
857 {
859 {
860 s11_values[j] += Scaling;
861 found = true;
862 }
863 }
864 }
865
867 const bool want_block = true;
868 for (unsigned b_i = 0; b_i < 2; b_i++)
869 {
870 for (unsigned b_j = 0; b_j < 2; b_j++)
871 {
872 required_blocks[b_i][b_j].select_block(b_i, b_j, want_block);
873 }
874 }
875
876 required_blocks[1][1].set_replacement_block_pt(s11_pt);
877
878 CRDoubleMatrix s_prec_pt = this->get_concatenated_block(required_blocks);
879
880 delete s11_pt;
881 s11_pt = 0;
882
883 // setup the preconditioner
885 {
886 Preconditioner_pt = (*Subsidiary_preconditioner_function_pt)();
887 }
888 else
889 {
892 }
894 }
895
896 //=============================================================================
897 /// Apply the preconditioner.
898 //=============================================================================
908
909
910 ///////////////////////////////////////////////////////////////////////////////
911 ///////////////////////////////////////////////////////////////////////////////
912 ///////////////////////////////////////////////////////////////////////////////
913
914
915 //=============================================================================
916 /// clean up the memory
917 //=============================================================================
920 {
921 // number of block types
923
924 // delete diagonal blocks
925 for (unsigned i = 0; i < n_block; i++)
926 {
929 if (Method == 1)
930 {
931 for (unsigned j = i + 1; j < n_block; j++)
932 {
935 }
936 }
937 else if (Method == 2)
938 {
939 for (unsigned j = 0; j < i; j++)
940 {
943 }
944 }
945 }
946
947 // clean up the block preconditioner
949 }
950
951 //=============================================================================
952 /// Setup the preconditioner.
953 //=============================================================================
955 {
956 // clean the memory
957 this->clean_up_memory();
958
959 // determine the number of DOF types
960 unsigned n_dof_types = this->ndof_types();
961
962#ifdef PARANOID
963 // must be Dim*2 dof types
964 if (n_dof_types % 2 != 0)
965 {
966 std::ostringstream error_message;
967 error_message << "This preconditioner requires DIM*3 types of DOF";
968 throw OomphLibError(
970 }
971#endif
972
973 // store the dimension of the problem
974 unsigned dim = n_dof_types / 2;
975
976 // assemble the dof to block lookup scheme
978 for (unsigned d = 0; d < dim; d++)
979 {
980 dof_to_block_map[d] = d;
981 dof_to_block_map[d + dim] = d;
982 }
983
984 // setup the blocks look up schemes
985 this->block_setup(dof_to_block_map);
986
987 // Storage for the diagonal block preconditioners
989
990 // storage for the off diagonal matrix vector products
991 Off_diagonal_matrix_vector_products.resize(dim, dim, 0);
992
993 // setup the subsidiary preconditioners
994 for (unsigned d = 0; d < dim; d++)
995 {
997 dof_list[0] = d;
998 dof_list[1] = d + dim;
999
1003 ->turn_into_subsidiary_block_preconditioner(this, dof_list);
1005 {
1007 ->set_subsidiary_preconditioner_function(
1009 }
1011
1012 Diagonal_block_preconditioner_pt[d]->Preconditioner::setup(matrix_pt());
1013
1014 // the preconditioning method.\n
1015 // 0 - block diagonal\n
1016 // 1 - upper triangular\n
1017 // 2 - lower triangular\n
1018 // next setup the off diagonal mat vec operators if required
1019 if (Method == 1 || Method == 2)
1020 {
1021 unsigned l = d + 1;
1022 unsigned u = dim;
1023 if (Method == 2)
1024 {
1025 l = 0;
1026 u = d;
1027 }
1028 for (unsigned j = l; j < u; j++)
1029 {
1031 this->get_block(d, j, *block_matrix_pt);
1033 // Off_diagonal_matrix_vector_products(d,j)->setup(block_matrix_pt);
1036
1037 delete block_matrix_pt;
1038 block_matrix_pt = 0;
1039 }
1040 }
1041 } // setup the subsidiary preconditioner.
1042 } // preconditioner setup
1043
1044 //=============================================================================
1045 /// Apply preconditioner to r
1046 //=============================================================================
1049 {
1050 // copy r
1052
1053 unsigned n_block;
1054
1055 // Cache umber of block types (also the spatial DIM)
1056 n_block = this->nblock_types();
1057
1058 // loop parameters
1059 int start = n_block - 1;
1060 int end = -1;
1061 int step = -1;
1062 if (Method != 1)
1063 {
1064 start = 0;
1065 end = n_block;
1066 step = 1;
1067 }
1068
1069 // the preconditioning method.
1070 // 0 - block diagonal
1071 // 1 - upper triangular
1072 // 2 - lower triangular
1073 //
1074 // loop over the DIM
1075 //
1076 // For Method = 0 or 2 (diagonal, lower)
1077 // start = 2, end = -1, step = -1
1078 // i = 2,1,0
1079 //
1080 // For Method = 1 (upper)
1081 // start = 0, end = 3 step = 1
1082 // i = 0, 1, 2
1083 for (int i = start; i != end; i += step)
1084 {
1085 // solve
1086 Diagonal_block_preconditioner_pt[i]->preconditioner_solve(r, z);
1087
1088 // if upper or lower triangular
1089 if (Method != 0)
1090 {
1091 // substitute
1092 //
1093 for (int j = i + step; j != end; j += step)
1094 {
1095 DoubleVector x;
1096 this->get_block_vector(i, z, x);
1097 DoubleVector y;
1098 Off_diagonal_matrix_vector_products(j, i)->multiply(x, y);
1099 x.clear();
1100 this->get_block_vector(j, r, x);
1101 x -= y;
1102 this->return_block_vector(j, x, r);
1103 } // substitute
1104 } // if upper or lower
1105 } // for loop over DIM
1106 } // Block preconditioner solve
1107} // namespace oomph
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...
unsigned ndof_types_in_mesh(const unsigned &i) const
Return the number of DOF types in mesh i. WARNING: This should only be used by the upper-most master ...
CRDoubleMatrix get_concatenated_block(const VectorMatrix< BlockSelector > &selected_block)
Returns a concatenation of the block matrices specified by the argument selected_block....
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 clear_block_preconditioner_base()
Clears all BlockPreconditioner data. Called by the destructor and the block_setup(....
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...
bool is_master_block_preconditioner() const
Return true if this preconditioner is the master block preconditioner.
unsigned nblock_types() const
Return the number of block types.
void set_replacement_dof_block(const unsigned &block_i, const unsigned &block_j, CRDoubleMatrix *replacement_dof_block_pt)
Set replacement dof-level blocks. Only dof-level blocks can be set. This is important due to how the ...
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 ...
unsigned ndof_types() const
Return the total number of DOF types.
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,...
void set_nmesh(const unsigned &n)
Specify the number of meshes required by this block preconditioner. Note: elements in different meshe...
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...
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
unsigned nrow_local() const
access function for the num of local rows on this processor.
A vector in the mathematical sense, initially developed for linear algebra type applications....
double * values_pt()
access function to the underlying values
void clear()
wipes the DoubleVector
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
An Preconditioner class using the suite of Hypre preconditioners to allow.
A preconditioner for performing inner iteration preconditioner solves. The template argument SOLVER s...
static OomphCommunicator * communicator_pt()
access to global communicator. This is the oomph-lib equivalent of MPI_COMM_WORLD
Matrix-based diagonal preconditioner.
Matrix vector product helper class - primarily a wrapper to Trilinos's Epetra matrix vector product m...
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...
virtual const OomphCommunicator * comm_pt() const
Get function for comm pointer.
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...
void lagrange_multiplier_preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the lagrange multiplier subsidiary preconditioner.
Preconditioner * Elastic_preconditioner_pt
storage for the preconditioner for the solid system
Vector< Preconditioner * > Lagrange_multiplier_preconditioner_pt
lagrange multiplier preconditioner pt
bool Use_inf_norm_of_s_scaling
boolean indicating whether the inf-norm of S should be used as scaling. Default = true;
SubsidiaryPreconditionerFctPt Elastic_subsidiary_preconditioner_function_pt
The solid subsidiary preconditioner function pointer.
Elastic_preconditioner_type E_preconditioner_type
An unsigned indicating which method should be used for preconditioning the solid component.
void elastic_preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the elastic subsidiary preconditioner.
double Scaling
The scaling. Defaults to infinity norm of S.
SubsidiaryPreconditionerFctPt Lagrange_multiplier_subsidiary_preconditioner_function_pt
The Lagrange multiplier subsidary preconditioner function pointer.
Mesh * Elastic_mesh_pt
Pointer to the mesh containing the solid elements.
Mesh * Lagrange_multiplier_mesh_pt
Pointer to the mesh containing the Lagrange multiplier elements.
A helper class for PseudoElasticPreconditioner. Note that this is NOT actually a functioning precondi...
Subsidiary helper preconditioner for the PseudoElasticPreconditioner. Required for block precondition...
void preconditioner_solve(const DoubleVector &res, DoubleVector &z)
Apply preconditioner to r.
SubsidiaryPreconditionerFctPt Subsidiary_preconditioner_function_pt
The SubisidaryPreconditionerFctPt.
Vector< PseudoElasticPreconditionerSubsidiaryPreconditionerOld * > Diagonal_block_preconditioner_pt
Vector of SuperLU preconditioner pointers for storing the preconditioners for each diagonal block.
unsigned Method
the preconditioning method. 0 - block diagonal 1 - upper triangular 2 - lower triangular
DenseMatrix< MatrixVectorProduct * > Off_diagonal_matrix_vector_products
Matrix of matrix vector product operators for the off diagonals.
Subsidiary helper preconditioner for the PseudoElasticPreconditioner. Required to construct the augme...
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the preconditioner.
SubsidiaryPreconditionerFctPt Subsidiary_preconditioner_function_pt
the SubisidaryPreconditionerFctPt
SubsidiaryPreconditionerFctPt Lagrange_multiplier_subsidiary_preconditioner_function_pt
The Lagrange multiplier subsidiary preconditioner function pointer.
void lagrange_multiplier_preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the lagrange multiplier subsidiary preconditioner.
bool Use_inf_norm_of_s_scaling
boolean indicating whether the inf-norm of S should be used as scaling. Default = true;
double Scaling
The scaling. Defaults to infinity norm of S.
Vector< Preconditioner * > Lagrange_multiplier_preconditioner_pt
lagrange multiplier preconditioner pt
void elastic_preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply the elastic subsidiary preconditioner.
Mesh * Lagrange_multiplier_mesh_pt
Pointer to the mesh containing the Lagrange multiplier elements.
SubsidiaryPreconditionerFctPt Elastic_subsidiary_preconditioner_function_pt
The solid subsidiary preconditioner function pointer.
Mesh * Elastic_mesh_pt
Pointer to the mesh containing the solid elements.
Elastic_preconditioner_type E_preconditioner_type
An unsigned indicating which method should be used for preconditioning the solid component.
Preconditioner * Elastic_preconditioner_pt
storage for the preconditioner for the solid system
unsigned Dim
the dimension of the problem
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
An interface to the Trilinos AztecOO classes allowing it to be used as an Oomph-lib LinearSolver....
An interface to the Trilinos ML class - provides a function to construct a serial ML object,...
double inf_norm(const DenseMatrix< CRDoubleMatrix * > &matrix_pt)
Compute infinity (maximum) norm of sub blocks as if it was one matrix.
Definition matrices.cc:3731
Preconditioner * create_exact_preconditioner()
Factory function to create suitable exact preconditioner.
Preconditioner * get_lagrange_multiplier_preconditioner()
CG with diagonal preconditioner for the lagrange multiplier subsidiary linear systems.
Preconditioner * get_elastic_preconditioner()
AMG w/ GS smoothing for the augmented elastic subsidiary linear systems – calls Hypre version to stay...
Preconditioner * get_elastic_preconditioner_trilinos_ml()
TrilinosML smoothing for the augmented elastic subsidiary linear systems.
Preconditioner * get_elastic_preconditioner_hypre()
AMG w/ GS smoothing for the augmented elastic subsidiary linear systems.
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).