solid_preconditioners.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//====================================================================
27
28
29namespace oomph
30{
31 //===========================================================================
32 /// Setup the least-squares commutator solid preconditioner. This
33 /// extracts blocks corresponding to the position/displacement and pressure
34 /// unknowns, creates the matrices actually needed in the application of the
35 /// preconditioner and deletes what can be deleted... Note that
36 /// this preconditioner needs a CRDoubleMatrix.
37 //============================================================================
39 {
40 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
41 // NOTE: In the interest of minimising memory usage, several containers
42 // are recycled, therefore their content/meaning changes
43 // throughout this function. The code is carefully annotated
44 // but you'll have to read it line by line!
45 //%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
46
47 // make sure any old data is deleted
49
50#ifdef PARANOID
51 // paranoid check that the solid mesh pt has been set
52 if (Solid_mesh_pt == 0)
53 {
54 std::ostringstream error_message;
55 error_message << "The solid mesh pointer must be set.\n"
56 << "Use method set_solid_mesh(...)";
57 throw OomphLibError(
59 }
60#endif
61
62 // set the mesh
63 this->set_mesh(0, Solid_mesh_pt);
64
65 // Get blocks
66 // ----------
67
68 // Set up block look up schemes (done automatically in the
69 // BlockPreconditioner base class, based on the information
70 // provided in the block-preconditionable elements in the problem)
71
72 // this preconditioner has two types of block:
73 // type 0: displacement/positions - corresponding to DOFs 0 to n-2
74 // type 1: pressure - corresponding to DOF n-1
76 unsigned ndof_types = 0;
78 {
79 ndof_types = this->ndof_types();
80 }
81 else
82 {
83 ndof_types = this->ndof_types_in_mesh(0);
84 }
87 this->block_setup(dof_to_block_map);
90 if (Doc_time)
91 {
92 oomph_info << "Time for block_setup(...) [sec]: " << block_setup_time
93 << "\n";
94 }
95
96 // determine whether the F preconditioner is a block preconditioner (and
97 // therefore a subsidiary preconditioner)
102 {
104 }
105
106 // Get B (the divergence block)
109 this->get_block(1, 0, *b_pt);
111 if (Doc_time)
112 {
114 oomph_info << "Time to get B [sec]: " << get_B_time << "\n";
115 }
116
117 // the pointer for f
119
120 // the pointer for the p matrix
122
123 // the pointer for bt
125
126 // if BFBt is to be formed
128 {
129 // If using scaling then replace B with Bt
132 {
133 // Assemble the ivmm diagonal
137 if (Doc_time)
138 {
139 double ivmm_assembly_time =
141 oomph_info << "Time to assemble inverse mass matrix [sec]: "
142 << ivmm_assembly_time << "\n";
143 }
144
145 // assemble BQ (stored in B)
148 b_pt->multiply(*ivmm_pt, *temp_matrix_pt);
149 delete b_pt;
150 b_pt = 0;
153 if (Doc_time)
154 {
156 oomph_info << "Time to generate BQ [sec]: " << t_BQ_time << std::endl;
157 }
158 }
159
160 // Get Bt
162 this->get_block(0, 1, *bt_pt);
164 if (Doc_time)
165 {
167 oomph_info << "Time to get Bt [sec]: " << t_get_Bt_time << std::endl;
168 }
169
170 // now form the P matrix by multiplying B (which if using scaling will be
171 // BQ) with Bt
173 b_pt->multiply(*bt_pt, *p_matrix_pt);
175 if (Doc_time)
176 {
177 double t_P_time = t_P_finish - t_P_start;
178 oomph_info << "Time to generate P matrix [sec]: " << t_P_time
179 << std::endl;
180 }
181
182 // Multiply auxiliary matrix by diagonal of mass matrix if
183 // required
185 {
188 ivmm_pt->multiply(*bt_pt, *temp_matrix_pt);
189 delete bt_pt;
190 bt_pt = 0;
193
194 // Output times
195 if (Doc_time)
196 {
198 oomph_info << "Time to generate QBt [sec]: " << t_QBt_time
199 << std::endl;
200 }
201 }
202
203 // Clean up memory
204 delete ivmm_pt;
205
206 // get block 0 0
208 this->get_block(0, 0, *f_pt);
210 if (Doc_time)
211 {
213 oomph_info << "Time to get F [sec]: " << t_get_F_time << std::endl;
214 }
215
216 // Auxiliary matrix for intermediate results
219 f_pt->multiply(*bt_pt, *aux_matrix_pt);
221 if (Doc_time)
222 {
224 oomph_info << "Time to generate FQBt [sec]: " << t_aux_time
225 << std::endl;
226 }
227
228 // can now delete Bt (or QBt if using scaling)
229 delete bt_pt;
230 bt_pt = 0;
231
232 // and if F requires a block preconditioner then we can delete F
234 {
235 delete f_pt;
236 }
237
238 // now form BFBt
241 b_pt->multiply(*aux_matrix_pt, *e_matrix_pt);
242 delete aux_matrix_pt;
243 delete b_pt;
245 if (Doc_time)
246 {
248 oomph_info << "Time to generate E (B*(F*Bt)) [sec]: " << t_E_time
249 << std::endl;
250 }
253 // E_mat_vec_pt->setup(e_matrix_pt);
254 this->setup_matrix_vector_product(E_mat_vec_pt, e_matrix_pt, 1);
256 delete e_matrix_pt;
257 if (Doc_time)
258 {
260 oomph_info << "Time to build E (BFBt) matrix vector operator E [sec]: "
261 << t_E_time << std::endl;
262 }
263
264 // and rebuild Bt
266 bt_pt = new CRDoubleMatrix;
267 this->get_block(0, 1, *bt_pt);
269 if (Doc_time)
270 {
272 oomph_info << "Time to get Bt [sec]: " << t_get_Bt_time << std::endl;
273 }
274 }
275
276
277 /////////////////////////////////////////////////////////////////////////////
278
279
280 // otherwise we are not forming BFBt
281 else
282 {
283 // get the inverse mass matrix (Q)
286 {
290 if (Doc_time)
291 {
292 double ivmm_assembly_time =
294 oomph_info << "Time to assemble Q (inverse diagonal "
295 << "mass matrix) [sec]: " << ivmm_assembly_time << "\n";
296 }
297 }
298
299 // Get Bt
301 this->get_block(0, 1, *bt_pt);
303 if (Doc_time)
304 {
306 oomph_info << "Time to get Bt [sec]: " << t_get_Bt_time << std::endl;
307 }
308
309 // next QBt
311 {
314 ivmm_pt->multiply(*bt_pt, *qbt_pt);
315 delete bt_pt;
316 bt_pt = 0;
317 bt_pt = qbt_pt;
319 if (Doc_time)
320 {
322 oomph_info << "Time to generate QBt [sec]: " << t_QBt_time
323 << std::endl;
324 }
325 delete ivmm_pt;
326 }
327
328 // form P
330 b_pt->multiply(*bt_pt, *p_matrix_pt);
332 if (Doc_time)
333 {
335 oomph_info << "Time to generate P [sec]: " << t_p_time << std::endl;
336 }
337 delete b_pt;
338 b_pt = 0;
339
340 // build the matvec operator for QBt
343 // QBt_mat_vec_pt->setup(bt_pt);
346 if (Doc_time)
347 {
349 oomph_info << "Time to build QBt matrix vector operator [sec]: "
350 << t_p_time << std::endl;
351 }
352 delete bt_pt;
353 bt_pt = 0;
354
355 // get F
357 this->get_block(0, 0, *f_pt);
359 if (Doc_time)
360 {
362 oomph_info << "Time to get F [sec]: " << t_get_F_time << std::endl;
363 }
364
365 // form the matrix vector product helper
368 // F_mat_vec_pt->setup(f_pt);
371 if (Doc_time)
372 {
374 oomph_info << "Time to build F Matrix Vector Operator [sec]: "
375 << t_F_MV_time << std::endl;
376 }
377
378 // if F is a block preconditioner then we can delete the F matrix
380 {
381 delete f_pt;
382 f_pt = 0;
383 }
384
385 // and rebuild Bt
387 bt_pt = new CRDoubleMatrix;
388 this->get_block(0, 1, *bt_pt);
390 if (Doc_time)
391 {
393 oomph_info << "Time to get Bt [sec]: " << t_get_Bt_time << std::endl;
394 }
395 }
396
397
398 /////////////////////////////////////////////////////////////////////////////
399
400
401 // form the matrix vector operator for Bt
404 // Bt_mat_vec_pt->setup(bt_pt);
407 if (Doc_time)
408 {
410 oomph_info << "Time to build Bt Matrix Vector Operator [sec]: "
411 << t_Bt_MV_time << std::endl;
412 }
413 delete bt_pt;
414 bt_pt = 0;
415
416 // if the P preconditioner has not been setup
417 if (P_preconditioner_pt == 0)
418 {
422 }
423
424 // std::stringstream junk;
425 // junk << "p_matrix" << comm_pt()->my_rank()
426 // << ".dat";
427 // p_matrix_pt->sparse_indexed_output_with_offset(junk.str());
428 // oomph_info << "Done output of " << junk.str() << std::endl;
429
430 // Setup the preconditioner for the Pressure matrix
433 delete p_matrix_pt;
434 p_matrix_pt = 0;
435 p_matrix_pt = 0;
437 if (Doc_time)
438 {
440 oomph_info << "P sub-preconditioner setup time [sec]: " << t_p_prec_time
441 << "\n";
442 }
443
444 // Set up solver for solution of system with momentum matrix
445 // ----------------------------------------------------------
446
447 // if the F preconditioner has not been setup
448 if (F_preconditioner_pt == 0)
449 {
453 }
454
455 // if F is a block preconditioner
458 {
459 unsigned ndof_types = this->ndof_types();
460 ndof_types--;
462 for (unsigned i = 0; i < ndof_types; i++)
463 {
464 dof_map[i] = i;
465 }
466 F_block_preconditioner_pt->turn_into_subsidiary_block_preconditioner(
467 this, dof_map);
469 }
470 // otherwise F is not a block preconditioner
471 else
472 {
474 delete f_pt;
475 f_pt = 0;
476 }
478 if (Doc_time)
479 {
481 oomph_info << "F sub-preconditioner setup time [sec]: " << t_f_prec_time
482 << "\n";
483 }
484
485 // Remember that the preconditioner has been setup so
486 // the stored information can be wiped when we
487 // come here next...
489 }
490
491
492 //=======================================================================
493 /// Apply preconditioner to r.
494 //=======================================================================
496 const DoubleVector& r, DoubleVector& z)
497 {
498#ifdef PARANOID
500 {
501 std::ostringstream error_message;
502 error_message << "setup must be called before using preconditioner_solve";
503 throw OomphLibError(
505 }
506 if (z.built())
507 {
508 if (z.nrow() != r.nrow())
509 {
510 std::ostringstream error_message;
511 error_message << "The vectors z and r must have the same number of "
512 << "of global rows";
513 throw OomphLibError(error_message.str(),
516 }
517 }
518#endif
519
521 double t_start = TimingHelpers::timer();
522 double t_end = 0;
523
524 // if z is not setup then give it the same distribution
525 if (!z.built())
526 {
527 z.build(r.distribution_pt(), 0.0);
528 }
529
530 // Step 1 - apply approximate Schur inverse to pressure unknowns (block 1)
531 // -----------------------------------------------------------------------
532
533 // Working vectors
536
537 // Copy pressure values from residual vector to temp_vec:
538 // Loop over all entries in the global vector (this one
539 // includes displacement/position and pressure dofs in some random fashion)
540 this->get_block_vector(1, r, temp_vec);
541
542
543 if (Doc_time)
544 {
546 oomph_info << "LSC prec solve: Time for get block vector: "
547 << t_end - t_start << std::endl;
549 }
550
551 // NOTE: The vector temp_vec now contains the vector r_p.
552
553 // Solve first pressure Poisson system
554#ifdef PARANOID
555 // check a solver has been set
556 if (P_preconditioner_pt == 0)
557 {
558 std::ostringstream error_message;
559 error_message << "P_preconditioner_pt has not been set.";
560 throw OomphLibError(
562 }
563#endif
564
565 // use some Preconditioner's preconditioner_solve function
567
568
569 if (Doc_time)
570 {
572 oomph_info << "LSC prec solve: First P solve [nrow="
573 << P_preconditioner_pt->nrow() << "]: " << t_end - t_start
574 << std::endl;
576 }
577
578
579 // NOTE: The vector another_temp_vec now contains the vector P^{-1} r_p
580
581 // Multiply another_temp_vec by matrix E and stick the result into temp_vec
582 temp_vec.clear();
584 {
586 }
587 else
588 {
590 another_temp_vec.clear();
592 temp_vec.clear();
594 }
595
596
597 if (Doc_time)
598 {
600 oomph_info << "LSC prec solve: E matrix vector product: "
601 << t_end - t_start << std::endl;
603 }
604
605 // NOTE: The vector temp_vec now contains E P^{-1} r_p
606
607 // Solve second pressure Poisson system using preconditioner_solve
608 another_temp_vec.clear();
610
611
612 if (Doc_time)
613 {
615 oomph_info << "LSC prec solve: Second P solve [nrow="
616 << P_preconditioner_pt->nrow() << "]: " << t_end - t_start
617 << std::endl;
619 }
620
621
622 // NOTE: The vector another_temp_vec now contains z_p = P^{-1} E P^{-1} r_p
623 // as required (apart from the sign which we'll fix in the
624 // next step.
625
626 // Now copy another_temp_vec (i.e. z_p) back into the global vector z.
627 // Loop over all entries in the global results vector z:
628 temp_vec.build(another_temp_vec.distribution_pt(), 0.0);
631
632 // Step 2 - apply preconditioner to displacement/positon unknowns (block 0)
633 // ------------------------------------------------------------------------
634
635 // Recall that another_temp_vec (computed above) contains the
636 // negative of the solution of the Schur complement systen, -z_p.
637 // Multiply by G (stored in Block_matrix_pt(0,1) and store
638 // result in temp_vec (vector resizes itself).
639 temp_vec.clear();
641
642
643 if (Doc_time)
644 {
646 oomph_info << "LSC prec solve: G matrix vector product: "
647 << t_end - t_start << std::endl;
649 }
650
651 // NOTE: temp_vec now contains -G z_p
652
653 // The vector another_temp_vec is no longer needed -- re-use it to store
654 // displacement/position quantities:
655 another_temp_vec.clear();
656
657 // Loop over all enries in the global vector and find the
658 // entries associated with the displacement/position:
661
662 // NOTE: The vector another_temp_vec now contains r_u - G z_p
663
664 // Solve momentum system
665#ifdef PARANOID
666 // check a solver has been set
667 if (F_preconditioner_pt == 0)
668 {
669 std::ostringstream error_message;
670 error_message << "F_preconditioner_pt has not been set.";
671 throw OomphLibError(
673 }
674#endif
675
676 // use some Preconditioner's preconditioner solve
677 // and return
679 {
682 }
683 else
684 {
687 }
688
689 if (Doc_time)
690 {
692 oomph_info << "LSC prec solve: F solve [nrow="
693 << P_preconditioner_pt->nrow() << "]: " << t_end - t_start
694 << std::endl;
695 oomph_info << "LSC prec solve: Overall " << t_end - t_start_overall
696 << std::endl;
697 }
698 }
699
700
701 //========================================================================
702 /// Helper function to assemble the diagonal of the
703 /// mass matrix from the elemental contributions defined in
704 /// SolidElementWithDiagonalMassMatrix::get_mass_matrix_diagonal(...).
705 //========================================================================
708 {
709 // determine the rows required by this processor
710 unsigned first_row = this->block_distribution_pt(0)->first_row();
711 unsigned nrow_local = this->block_distribution_pt(0)->nrow_local();
712 unsigned nrow = this->block_distribution_pt(0)->nrow();
713
714 // create storage for the diagonals
715 double* m_values = new double[nrow_local];
716 for (unsigned i = 0; i < nrow_local; i++)
717 {
718 m_values[i] = 0;
719 }
720
721 // if the problem is distributed
722 bool distributed = false;
723#ifdef OOMPH_HAS_MPI
725 {
726 distributed = true;
727 }
728#endif
729
730 // next we get the diagonal mass matrix data
731 if (distributed)
732 {
733#ifdef OOMPH_HAS_MPI
734 // the number of processors
735 unsigned nproc = this->comm_pt()->nproc();
736
737 // and my rank
738 unsigned my_rank = this->comm_pt()->my_rank();
739
740 // determine the rows for which we have lookup rows
741 // if the problem is NOT distributed then we only classify global equation
742 // on this processor to avoid duplication (as every processor holds
743 // every element)
744 unsigned first_lookup_row = 0;
745 unsigned last_lookup_row = 0;
746 first_lookup_row = this->master_distribution_pt()->first_row();
747 last_lookup_row =
748 first_lookup_row + this->master_distribution_pt()->nrow_local() - 1;
749
750 // find number of local elements
751 unsigned n_el = Solid_mesh_pt->nelement();
752
753 // the diagonal mass matrix contributions that have been
754 // classified and should be sent to another processor
756
757 // the corresponding block indices
759
760 // the maitrix contributions that cannot be classified by this processor
761 // and therefore must be sent to another for classification
764
765 // the corresponding global indices that require classification
767
768 // get the master distribution pt
771
772 // get the displacement/position distribution pt
774 this->block_distribution_pt(0);
775
776 // get the contribution for each element
777 for (unsigned e = 0; e < n_el; e++)
778 {
779 // check that the element is not halo d
781 {
782 // find number of degrees of freedom in the element
783 // (this is slightly too big because it includes the
784 // pressure dofs but this doesn't matter)
785 unsigned el_dof = Solid_mesh_pt->element_pt(e)->ndof();
786
787 // allocate local storage for the element's contribution to the
788 // mass matrix diagonal
790
794
795
796 if (cast_el_pt == 0)
797 {
798#ifdef PARANOID
799 std::ostringstream error_message;
800 error_message << "Failed cast to "
801 << "SolidElementWithDiagonalMassMatrix*\n"
802 << "Element is of type: "
803 << typeid(*(Solid_mesh_pt->element_pt(e))).name()
804 << "\n"
805 << typeid(Solid_mesh_pt->element_pt(e)).name()
806 << std::endl;
807 OomphLibWarning(error_message.str(),
808 "PressureBasedSolidLSCPreconditioner::assemble_"
809 "mass_matrix_diagonal()",
811#endif
812 }
813 else
814 {
815 cast_el_pt->get_mass_matrix_diagonal(el_vmm_diagonal);
816 }
817
818 // get the contribution for each dof
819 for (unsigned i = 0; i < el_dof; i++)
820 {
821 // Get the equation number
822 unsigned eqn_number = Solid_mesh_pt->element_pt(e)->eqn_number(i);
823
824 // if I have lookup information on this processor
825 if (eqn_number >= first_lookup_row && eqn_number <= last_lookup_row)
826 {
827 // bypass non displacement/position DOFs
828 if (this->block_number(eqn_number) == 0)
829 {
830 // get the index in the block
831 unsigned index = this->index_in_block(eqn_number);
832
833 // determine which processor requires the block index
834 for (unsigned p = 0; p < nproc; p++)
835 {
836 if (index >= displ_dist_pt->first_row(p) &&
837 (index < (displ_dist_pt->first_row(p) +
838 displ_dist_pt->nrow_local(p))))
839 {
840 // if it is required by this processor then add the
841 // contribution
842 if (p == my_rank)
843 {
844 m_values[index - first_row] += el_vmm_diagonal[i];
845 }
846 // other wise store it for communication
847 else
848 {
851 classified_indices_send[p].push_back(index);
852 }
853 }
854 }
855 }
856 }
857
858 // if we do not have the lookup information on this processor
859 // then we send the mass matrix contribution to a processor
860 // which we know does have the lookup information
861 // the assumption: the processor for which the master block
862 // preconditioner distribution will definitely hold the lookup
863 // data for eqn_number (although others may)
864 else if (any_mesh_distributed())
865 {
866 // determine which processor requires the block index
867 unsigned p = 0;
868 while (!(eqn_number >= master_distribution_pt->first_row(p) &&
869 (eqn_number < (master_distribution_pt->first_row(p) +
870 master_distribution_pt->nrow_local(p)))))
871 {
872 p++;
873 }
874
875 // store the data
877 unclassified_indices_send[p].push_back(eqn_number);
878 }
879 }
880 }
881 }
882
883 // next the unclassified contributions are communicated to
884 // processors that can classify them
885
886 // first determine how many unclassified rows are to be sent to
887 // each processor
888 unsigned* n_unclassified_send = new unsigned[nproc];
889 for (unsigned p = 0; p < nproc; p++)
890 {
891 if (p == my_rank)
892 {
894 }
895 else
896 {
898 }
899 }
900
901 // then all-to-all number of unclassified to be sent / recv
902 unsigned* n_unclassified_recv = new unsigned[nproc];
904 1,
907 1,
909 this->comm_pt()->mpi_comm());
910
911 // the base displacement for the sends
914
915 // allocate storage for the data to be recieved
916 // and post the sends and recvs
922 for (unsigned p = 0; p < nproc; p++)
923 {
924 if (p != my_rank)
925 {
926 // recv
927 if (n_unclassified_recv[p] > 0)
928 {
930 new double[n_unclassified_recv[p]];
932
933 // data for the struct data type
936 int recv_sz[2];
937
938 // contributions
945 recv_sz[0] = 1;
946
947 // indices
954 recv_sz[1] = 1;
955
956 // build the final recv type
961
962 // and recv
964 MPI_Irecv(
965 m_values, 1, final_recv_type, p, 0, comm_pt()->mpi_comm(), &req);
967 unclassified_recv_proc.push_back(p);
971 }
972
973 // send
974 if (n_unclassified_send[p] > 0)
975 {
976 // data for the struct data type
979 int send_sz[2];
980
981 // contributions
988 send_sz[0] = 1;
989
990 // indices
997 send_sz[1] = 1;
998
999 // build the final send type
1004
1005 // and send
1007 MPI_Isend(
1008 m_values, 1, final_send_type, p, 0, comm_pt()->mpi_comm(), &req);
1013 }
1014 }
1015 }
1016
1017 // next classify the data as it is received
1019 while (n_unclassified_recv_req > 0)
1020 {
1021 // get the processor number and remove the completed request
1022 // for the vector of requests
1023 int req_num;
1026 &req_num,
1028 unsigned p = unclassified_recv_proc[req_num];
1030 req_num);
1033
1034 // next classify the dofs
1035 // and store them for sending to other processors if required
1036 unsigned n_recv = n_unclassified_recv[p];
1037 for (unsigned i = 0; i < n_recv; i++)
1038 {
1039 unsigned eqn_number = unclassified_indices_recv[p][i];
1040 // bypass non displacement/position DOFs
1041 if (this->block_number(eqn_number) == 0)
1042 {
1043 // get the index in the block
1044 unsigned index = this->index_in_block(eqn_number);
1045
1046 // determine which processor requires the block index
1047 for (unsigned pp = 0; pp < nproc; pp++)
1048 {
1049 if (index >= displ_dist_pt->first_row(pp) &&
1050 (index < (displ_dist_pt->first_row(pp) +
1051 displ_dist_pt->nrow_local(pp))))
1052 {
1053 // if it is required by this processor then add the
1054 // contribution
1055 if (pp == my_rank)
1056 {
1057 m_values[index - first_row] +=
1059 }
1060 // other wise store it for communication
1061 else
1062 {
1065 classified_indices_send[pp].push_back(index);
1066 }
1067 }
1068 }
1069 }
1070 }
1071
1072 // clean up
1074 delete[] unclassified_indices_recv[p];
1075 }
1076 delete[] n_unclassified_recv;
1077
1078 // now all indices have been classified
1079
1080 // next the classified contributions are communicated to
1081 // processors that require them
1082
1083 // first determine how many classified rows are to be sent to
1084 // each processor
1085 unsigned* n_classified_send = new unsigned[nproc];
1086 for (unsigned p = 0; p < nproc; p++)
1087 {
1088 if (p == my_rank)
1089 {
1090 n_classified_send[p] = 0;
1091 }
1092 else
1093 {
1095 }
1096 }
1097
1098 // then all-to-all com number of classified to be sent / recv
1099 unsigned* n_classified_recv = new unsigned[nproc];
1101 1,
1104 1,
1106 this->comm_pt()->mpi_comm());
1107
1108 // allocate storage for the data to be recieved
1109 // and post the sends and recvs
1115 for (unsigned p = 0; p < nproc; p++)
1116 {
1117 if (p != my_rank)
1118 {
1119 // recv
1120 if (n_classified_recv[p] > 0)
1121 {
1123 classified_indices_recv[p] = new unsigned[n_classified_recv[p]];
1124
1125 // data for the struct data type
1128 int recv_sz[2];
1129
1130 // contributions
1135 &recv_displacements[0]);
1137 recv_sz[0] = 1;
1138
1139 // indices
1145 recv_sz[1] = 1;
1146
1147 // build the final recv type
1152
1153 // and recv
1155 MPI_Irecv(
1156 m_values, 1, final_recv_type, p, 0, comm_pt()->mpi_comm(), &req);
1157 classified_recv_requests.push_back(req);
1158 classified_recv_proc.push_back(p);
1162 }
1163
1164 // send
1165 if (n_classified_send[p] > 0)
1166 {
1167 // data for the struct data type
1170 int send_sz[2];
1171
1172 // contributions
1177 &send_displacements[0]);
1179 send_sz[0] = 1;
1180
1181 // indices
1186 &send_displacements[1]);
1188 send_sz[1] = 1;
1189
1190 // build the final send type
1195
1196 // and send
1198 MPI_Isend(
1199 m_values, 1, final_send_type, p, 0, comm_pt()->mpi_comm(), &req);
1200 classified_send_requests.push_back(req);
1204 }
1205 }
1206 }
1207
1208 // next classify the data as it is received
1210 while (n_classified_recv_req > 0)
1211 {
1212 // get the processor number and remove the completed request
1213 // for the vector of requests
1214 int req_num;
1217 &req_num,
1219 unsigned p = classified_recv_proc[req_num];
1221 req_num);
1224
1225 // next classify the dofs
1226 // and store them for sending to other processors if required
1227 unsigned n_recv = n_classified_recv[p];
1228 for (unsigned i = 0; i < n_recv; i++)
1229 {
1232 }
1233
1234 // clean up
1236 delete[] classified_indices_recv[p];
1237 }
1238
1239 // wait for the unclassified sends to complete
1242 {
1246 }
1249 delete[] n_unclassified_send;
1250
1251 // wait for the classified sends to complete
1253 if (n_classified_send_req > 0)
1254 {
1258 }
1259 delete[] classified_indices_send;
1261 delete[] n_classified_recv;
1262 delete[] n_classified_send;
1263#endif
1264 }
1265
1266 // or if the problem is not distributed
1267 else
1268 {
1269 // find number of elements
1270 unsigned n_el = Solid_mesh_pt->nelement();
1271
1272 // get the contribution for each element
1273 for (unsigned e = 0; e < n_el; e++)
1274 {
1275 // find number of degrees of freedom in the element
1276 // (this is slightly too big because it includes the
1277 // pressure dofs but this doesn't matter)
1278 unsigned el_dof = Solid_mesh_pt->element_pt(e)->ndof();
1279
1280 // allocate local storage for the element's contribution to the
1281 // mass matrix diagonal
1283
1287
1288 if (cast_el_pt == 0)
1289 {
1290#ifdef PARANOID
1291 // #pragma clang diagnostic push
1292 // #pragma clang diagnostic ignored
1293 // "-Wpotentially-evaluated-expression"
1294 std::ostringstream error_message;
1295 error_message << "Failed cast to "
1296 << "SolidElementWithDiagonalMassMatrix*\n"
1297 << "Element is of type: "
1298 << typeid(*(Solid_mesh_pt->element_pt(e))).name()
1299 << "\n"
1300 << typeid(Solid_mesh_pt->element_pt(e)).name()
1301 << std::endl;
1302 OomphLibWarning(error_message.str(),
1303 "PressureBasedSolidLSCPreconditioner::assemble_mass_"
1304 "matrix_diagonal()",
1306// #pragma clang diagnostic pop
1307#endif
1308 }
1309 else
1310 {
1311 cast_el_pt->get_mass_matrix_diagonal(el_vmm_diagonal);
1312 }
1313
1314 // get the contribution for each dof
1315 for (unsigned i = 0; i < el_dof; i++)
1316 {
1317 // Get the equation number
1318 unsigned eqn_number = Solid_mesh_pt->element_pt(e)->eqn_number(i);
1319
1320 // bypass non displacement/position DOFs
1321 if (this->block_number(eqn_number) == 0)
1322 {
1323 // get the index in the block
1324 unsigned index = this->index_in_block(eqn_number);
1325
1326 // if it is required on this processor
1327 if (index >= first_row && index < first_row + nrow_local)
1328 {
1329 m_values[index - first_row] += el_vmm_diagonal[i];
1330 }
1331 }
1332 }
1333 }
1334 }
1335
1336 // create column index and row start
1337 int* m_column_index = new int[nrow_local];
1338 int* m_row_start = new int[nrow_local + 1];
1339 for (unsigned i = 0; i < nrow_local; i++)
1340 {
1341 m_values[i] = 1 / m_values[i];
1343 m_row_start[i] = i;
1344 }
1346
1347 // build the matrix
1349 m_pt->build_without_copy(
1351
1352 // return the matrix;
1353 return m_pt;
1354 }
1355
1356 //=========================================================================
1357 /// Helper function to delete preconditioner data.
1358 //=========================================================================
1360 {
1362 {
1363 // delete matvecs
1364 delete Bt_mat_vec_pt;
1365 Bt_mat_vec_pt = 0;
1366
1367 if (!Form_BFBt_product)
1368 {
1369 delete F_mat_vec_pt;
1370 F_mat_vec_pt = 0;
1371 delete QBt_mat_vec_pt;
1372 QBt_mat_vec_pt = 0;
1373 }
1374 else
1375 {
1376 delete E_mat_vec_pt;
1377 E_mat_vec_pt = 0;
1378 }
1379
1380 // delete stuff from displacement/position solve
1382 {
1383 delete F_preconditioner_pt;
1385 }
1386
1387 // delete stuff from Schur complement approx
1389 {
1390 delete P_preconditioner_pt;
1392 }
1393 }
1394 }
1395} // namespace oomph
e
Definition cfortran.h:571
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 ...
const LinearAlgebraDistribution * master_distribution_pt() const
Access function to the distribution of the master preconditioner. If this preconditioner does not hav...
bool any_mesh_distributed() const
Check if any of the meshes are distributed. This is equivalent to problem.distributed() and is used a...
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...
int index_in_block(const unsigned &i_dof) const
Given a global dof number, returns the index in the block it belongs to. This is the overall index,...
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...
int block_number(const unsigned &i_dof) const
Return the block number corresponding to a global index i_dof.
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 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*,...
const LinearAlgebraDistribution * block_distribution_pt(const unsigned &b) const
Access function to the block distributions (const version).
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...
A class for compressed row matrices. This is a distributable object.
Definition matrices.h:888
void build_without_copy(const unsigned &ncol, const unsigned &nnz, double *value, int *column_index, int *row_start)
keeps the existing distribution and just matrix that is stored without copying the matrix data
Definition matrices.cc:1710
bool distributed() const
distribution is serial or distributed
unsigned nrow() const
access function to the number of global rows.
unsigned nrow_local() const
access function for the num of local rows on this processor.
unsigned first_row() const
access function for the first row on this processor
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 size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
bool is_halo() const
Is this element a halo?
Definition elements.h:1150
unsigned ndof() const
Return the number of equations/dofs in the element.
Definition elements.h:822
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Definition elements.h:691
Describes the distribution of a distributable linear algebra type object. Typically this is a contain...
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.
GeneralisedElement *& element_pt(const unsigned long &e)
Return pointer to element e.
Definition mesh.h:452
unsigned long nelement() const
Return number of elements in the mesh.
Definition mesh.h:598
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
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...
MatrixVectorProduct * QBt_mat_vec_pt
MatrixVectorProduct operator for QBt if BFBt is not to be formed.
bool Form_BFBt_product
indicates whether BFBt should be formed or the component matrices should be retained....
CRDoubleMatrix * assemble_mass_matrix_diagonal()
Helper function to assemble the diagonal of the mass matrix from the elemental contributions defined ...
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.
MatrixVectorProduct * E_mat_vec_pt
MatrixVectorProduct operator for E (BFBt) if BFBt is to be formed.
bool Doc_time
Set Doc_time to true for outputting results of timings.
Mesh * Solid_mesh_pt
the pointer to the mesh of block preconditionable solid elements.
bool F_preconditioner_is_block_preconditioner
Boolean indicating whether the momentum system preconditioner is a block preconditioner.
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to Vector r.
MatrixVectorProduct * Bt_mat_vec_pt
MatrixVectorProduct operator for Bt;.
bool Using_default_p_preconditioner
flag indicating whether the default P preconditioner is used
MatrixVectorProduct * F_mat_vec_pt
MatrixVectorProduct operator for F if BFBt is not to be formed.
void clean_up_memory()
Helper function to delete preconditioner data.
bool Using_default_f_preconditioner
flag indicating whether the default F preconditioner is used
Preconditioner * P_preconditioner_pt
Pointer to the 'preconditioner' for the pressure matrix.
bool P_matrix_using_scaling
Control flag is true if mass matrix diagonal scaling is used in the Schur complement approximation.
Pure virtual base class for elements that can be used with PressureBasedSolidLSCPreconditioner.
Definition elements.h:5197
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
Preconditioner * create_exact_preconditioner()
Factory function to create suitable exact preconditioner.
double timer()
returns the time in seconds after some point in past
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).
OomphInfo oomph_info
Single (global) instantiation of the OomphInfo object – this is used throughout the library as a "rep...