general_purpose_space_time_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-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// Oomph-lib headers
28
29// Header file
31
32////////////////////////////////////////////////////////////////////////////
33////////////////////////////////////////////////////////////////////////////
34////////////////////////////////////////////////////////////////////////////
35
36namespace oomph
37{
38 //============================================================================
39 /// Set up the exact preconditioner
40 //============================================================================
41 template<typename MATRIX>
43 {
44 // Clean the memory
45 this->clean_up_memory();
46
47 // Subsidiary preconditioners don't really need the meshes
48 if (this->is_master_block_preconditioner())
49 {
50#ifdef PARANOID
51 if (this->gp_nmesh() == 0)
52 {
53 std::ostringstream err_msg;
54 err_msg << "There are no meshes set.\n"
55 << "Did you remember to call add_mesh(...)?";
56 throw OomphLibError(
58 }
59#endif
60
61 // Set all meshes if this is master block preconditioner
62 this->gp_preconditioner_set_all_meshes();
63 }
64
65 // If we're meant to build silently
66 if (this->Silent_preconditioner_setup == true)
67 {
68 // Store the output stream pointer
69 this->Stream_pt = oomph_info.stream_pt();
70
71 // Now set the oomph_info stream pointer to the null stream to
72 // disable all possible output
74 }
75
76 // Set up the block look up schemes
77 this->gp_preconditioner_block_setup();
78
79 // Number of block types
80 unsigned n_block_types = this->nblock_types();
81
82 // Fill in any null subsidiary preconditioners
83 this->fill_in_subsidiary_preconditioners(n_block_types);
84
85 // The total time for setting up the subsidiary preconditioners
86 double t_subsidiary_setup_total = 0.0;
87
88 // Stores the blocks we want to extract
90
91 // Boolean indicating if we want the block, stored for readability.
92 const bool want_block = true;
93
94 // Loop over the block rows
95 for (unsigned i = 0; i < n_block_types; i++)
96 {
97 // Loop over the block columns
98 for (unsigned j = 0; j < n_block_types; j++)
99 {
100 // Indicate that we want this block
101 required_blocks[i][j].select_block(i, j, want_block);
102 }
103 } // for (unsigned i=0;i<n_block_types;i++)
104
105 // Get the start time
107
108 // Get the concatenated blocks
109 CRDoubleMatrix full_matrix = this->get_concatenated_block(required_blocks);
110
111 // The total time for extracting all the blocks from the "global" matrix
113
114 // Get the start time
116
117 // It's a block preconditioner so pass it a pointer to the global matrix!
118 this->Subsidiary_preconditioner_pt[0]->setup(&full_matrix);
119
120 // Update the timing total
123
124 // Remember that the preconditioner has been set up
125 Preconditioner_has_been_setup = true;
126
127 // If we're meant to build silently, reassign the oomph stream pointer
128 if (this->Silent_preconditioner_setup == true)
129 {
130 // Store the output stream pointer
131 oomph_info.stream_pt() = this->Stream_pt;
132
133 // Reset our own stream pointer
134 this->Stream_pt = 0;
135 }
136
137 // Tell the user
138 oomph_info << "Total block extraction time [sec]: " << t_extraction_total
139 << "\nTotal subsidiary preconditioner setup time [sec]: "
140 << t_subsidiary_setup_total << std::endl;
141
142 } // End of setup
143
144 //==========================================================================
145 /// Preconditioner solve for the exact preconditioner
146 //==========================================================================
147 template<typename MATRIX>
149 const DoubleVector& r, DoubleVector& z)
150 {
151 // Vector of vectors for each section of residual vector
153
154 // Rearrange the vector r into the vector of block vectors block_r
155 this->get_block_vectors(r, block_r);
156
157 // Allocate space for the properly rearranged RHS vector
159
160 // Concatenate the DoubleVectors together
162
163 // Allocate space for the rearranged solution vector
165
166 // Solve the whole system exactly
167 this->Subsidiary_preconditioner_pt[0]->preconditioner_solve(rhs_reordered,
169
170 // Vector of vectors for each section of solution vector (copy the RHS
171 // block vector to get the sizing and distributions right)
173
174 // Split the solution vector into blocks
176
177 // Copy the solution from the block vector block_z back into z
178 this->return_block_vectors(block_z, z);
179 } // End of preconditioner_solve
180
181 //============================================================================
182 /// Set up the block triangular preconditioner
183 //============================================================================
184 template<typename MATRIX>
186 {
187 // Clean the memory
188 this->clean_up_memory();
189
190 // Subsidiary preconditioners don't really need the meshes
191 if (this->is_master_block_preconditioner())
192 {
193#ifdef PARANOID
194 if (this->gp_nmesh() == 0)
195 {
196 std::ostringstream err_msg;
197 err_msg << "There are no meshes set.\n"
198 << "Did you remember to call add_mesh(...)?";
199 throw OomphLibError(
201 }
202#endif
203
204 // Set all meshes if this is master block preconditioner
205 this->gp_preconditioner_set_all_meshes();
206 }
207
208 // If we're meant to build silently
209 if (this->Silent_preconditioner_setup == true)
210 {
211 // Store the output stream pointer
212 this->Stream_pt = oomph_info.stream_pt();
213
214 // Now set the oomph_info stream pointer to the null stream to
215 // disable all possible output
217 }
218
219 // Set up the block look up schemes
220 this->gp_preconditioner_block_setup();
221
222 // Number of block types
223 unsigned n_block_types = this->nblock_types();
224
225 // Storage for the off diagonal matrix vector products
226 Off_diagonal_matrix_vector_products.resize(n_block_types, n_block_types, 0);
227
228 // Fill in any null subsidiary preconditioners
229 this->fill_in_subsidiary_preconditioners(n_block_types);
230
231 // The total time for extracting all the blocks from the "global" matrix
232 double t_extraction_total = 0.0;
233
234 // The total time for setting up the subsidiary preconditioners
235 double t_subsidiary_setup_total = 0.0;
236
237 // The total time for setting up the matrix-vector products
238 double t_mvp_setup_total = 0.0;
239
240 // Build the preconditioners and matrix vector products
241 for (unsigned i = 0; i < n_block_types; i++)
242 {
243 // See if the i-th subsidiary preconditioner is a block preconditioner
244 if (dynamic_cast<BlockPreconditioner<CRDoubleMatrix>*>(
245 this->Subsidiary_preconditioner_pt[i]) != 0)
246 {
247 // Get the start time
249
250 // It's a block preconditioner so pass it a pointer to the global
251 // matrix!
252 this->Subsidiary_preconditioner_pt[i]->setup(this->matrix_pt());
253
254 // Get the end time
256
257 // Update the timing total
260 }
261 // It's just a preconditioner so pass it a deep copy of the block!
262 else
263 {
264 // Get the start time
266
267 // Grab the i-th diagonal block
268 CRDoubleMatrix block_matrix = this->get_block(i, i);
269
270 // Get the end time
272
273 // Update the timing total
275
276 // Get the start time
278
279 // Set up the i-th subsidiary preconditioner with this block
280 this->Subsidiary_preconditioner_pt[i]->setup(&block_matrix);
281
282 // Get the end time
284
285 // Update the timing total
288 }
289
290 // The index of the lower column bound
291 unsigned l;
292
293 // The index of the upper column bound
294 unsigned u;
295
296 // If we're using the upper triangular form of the preconditioner
297 if (Upper_triangular)
298 {
299 // The lower bound column bound
300 l = i + 1;
301
302 // If the block bandwidth has been provided
303 if (Block_bandwidth >= 0)
304 {
305 // The upper bound is the last nonzero block column
306 u = std::min(n_block_types, l + Block_bandwidth);
307 }
308 // Default case; loop over all the block columns
309 else
310 {
311 // The upper bound is the last column
312 u = n_block_types;
313 }
314 }
315 // If we're using the lower triangular form of the preconditioner
316 else
317 {
318 // The upper bound column bound
319 u = i;
320
321 // If the block bandwidth has been provided
322 if (Block_bandwidth >= 0)
323 {
324 // The lower bound is the last nonzero block column. The explicit
325 // casting is a bit over the top but it's better to be safe than
326 // sorry...
327 l = std::max(0, int(int(u) - Block_bandwidth));
328 }
329 // Default case; loop over all the block columns
330 else
331 {
332 // The lower bound is the first column
333 l = 0;
334 }
335 } // for (unsigned i=0;i<nblock_types;i++)
336
337 // Loop over the chosen block columns
338 for (unsigned j = l; j < u; j++)
339 {
340 // Get the start time
342
343 // Get the (i,j)-th block matrix
344 CRDoubleMatrix block_matrix = this->get_block(i, j);
345
346 // Get the end time
348
349 // Update the timing total
351
352 // Get the start time
354
355 // Copy the block into a "multiplier" class. If trilinos is being
356 // used this should also be faster than oomph-lib's multiphys.
357 Off_diagonal_matrix_vector_products(i, j) = new MatrixVectorProduct();
358
359 // Set the damn thing up
360 this->setup_matrix_vector_product(
361 Off_diagonal_matrix_vector_products(i, j), &block_matrix, j);
362
363 // Get the end time
365
366 // Update the timing total
368 }
369 } // for (unsigned i=0;i<n_block_types;i++)
370
371 // Remember that the preconditioner has been set up
372 Preconditioner_has_been_setup = true;
373
374 // If we're meant to build silently, reassign the oomph stream pointer
375 if (this->Silent_preconditioner_setup == true)
376 {
377 // Store the output stream pointer
378 oomph_info.stream_pt() = this->Stream_pt;
379
380 // Reset our own stream pointer
381 this->Stream_pt = 0;
382 }
383
384 // Tell the user
385 oomph_info << "Total block extraction time [sec]: " << t_extraction_total
386 << "\nTotal subsidiary preconditioner setup time [sec]: "
388 << "\nTotal matrix-vector product setup time [sec]: "
389 << t_mvp_setup_total << std::endl;
390
391 } // End of setup
392
393 //=============================================================================
394 /// Preconditioner solve for the block triangular preconditioner
395 //=============================================================================
396 template<typename MATRIX>
398 const DoubleVector& r, DoubleVector& z)
399 {
400 // Cache number of block types
401 unsigned n_block = this->nblock_types();
402
403 // The start index
404 int start;
405
406 // The end index
407 int end;
408
409 // The max. bandwidth of the block matrix
410 int bandwidth_limit;
411
412 // The steps to take (=+1/-1) which is dependent on whether we're using
413 // the lower or upper triangular form of the preconditioner
414 int step;
415
416 // If we're using this as an block upper triangular preconditioner
417 if (Upper_triangular)
418 {
419 start = n_block - 1;
420 end = -1;
421 step = -1;
422 }
423 else
424 {
425 start = 0;
426 end = n_block;
427 step = 1;
428 }
429
430 // Storage for the iteration count. The entries will be updated as we loop
431 // over the blocks. Here -1 is used to indicate that a block was not solved
432 // by GMRES so we ignore it
434
435 // Storage for the number of blocks solved with GMRES
436 unsigned n_block_solved_with_gmres = 0;
437
438 // Vector of vectors for each section of residual vector
440
441 // Rearrange the vector r into the vector of block vectors block_r
442 this->get_block_vectors(r, block_r);
443
444 // Vector of vectors for the solution block vectors
446
447 // Loop over the block rows
448 for (int i = start; i != end; i += step)
449 {
450 // If the i-th subsidiary preconditioner is a block preconditioner
451 if (dynamic_cast<BlockPreconditioner<CRDoubleMatrix>*>(
452 this->Subsidiary_preconditioner_pt[i]) != 0)
453 {
454 // This is pretty inefficient but there is no other way to do this with
455 // block sub preconditioners as far as I can tell because they demand
456 // to have the full r and z vectors...
457 DoubleVector block_z_with_size_of_full_z(r.distribution_pt());
458 DoubleVector r_updated(r.distribution_pt());
459
460 // Construct the new r vector (with the update given by back subs
461 // below).
462 this->return_block_vectors(block_r, r_updated);
463
464 // Solve (blocking is handled inside the block preconditioner).
465 this->Subsidiary_preconditioner_pt[i]->preconditioner_solve(
467
468 // Extract this block's z vector because we need it for the next step
469 // (and throw away block_z_with_size_of_full_z).
470 this->get_block_vector(i, block_z_with_size_of_full_z, block_z[i]);
471 }
472 // If the subsidiary preconditioner is just a regular preconditioner
473 else
474 {
475 // Solve on the block
476 this->Subsidiary_preconditioner_pt[i]->preconditioner_solve(block_r[i],
477 block_z[i]);
478 }
479
480 // See if the i-th subsidiary preconditioner the GMRES block
481 // preconditioner
482 if (dynamic_cast<GMRESBlockPreconditioner*>(
483 this->Subsidiary_preconditioner_pt[i]) != 0)
484 {
485 // Update the iteration count
487 this->Subsidiary_preconditioner_pt[i])
488 ->iterations();
489
490 // Increment the counter
492 }
493
494 // If the block bandwidth has been provided
495 if (Block_bandwidth >= 0)
496 {
497 // If we're using this as an block upper triangular preconditioner
498 if (Upper_triangular)
499 {
500 bandwidth_limit = std::max(end - i, step + step * Block_bandwidth);
501 }
502 else
503 {
504 bandwidth_limit = std::min(end - i, step + step * Block_bandwidth);
505 }
506 }
507 // Default case; loop over all the block rows
508 else
509 {
511 }
512
513 // Substitute (over all the rows)
514 for (int j = i + step; j != i + bandwidth_limit; j += step)
515 {
516 // Allocate space for the matrix-vector product (MVP)
518
519 // Calculate the MVP
520 Off_diagonal_matrix_vector_products(j, i)->multiply(block_z[i], temp);
521
522 // Now update the RHS vector
523 block_r[j] -= temp;
524 }
525 } // for (int i=start;i!=end;i+=step)
526
527 // Copy the solution from the block vector block_z back into z
528 this->return_block_vectors(block_z, z);
529
530 // Only compute the iteration count statistics if we have some data
532 {
533 // Storage for the average iteration count
534 double n_iter_mean = 0.0;
535
536 // Storage for the variance of the iteration count
537 double n_iter_var = 0.0;
538
539 // Loop over the entries of the iteration count vector
540 for (unsigned i = 0; i < n_block; i++)
541 {
542 // If we used GMRES on this block
543 if (dynamic_cast<GMRESBlockPreconditioner*>(
544 this->Subsidiary_preconditioner_pt[i]) != 0)
545 {
546 // Sanity check: make sure there's a non-negative iteration count
547 if (iteration_count[i] < 0)
548 {
549 // Create an output stream
550 std::ostringstream error_message_stream;
551
552 // Create an error message
554 << "Iteration count should be non-negative but "
555 << "you have an iteration\ncount of " << iteration_count[i] << "!"
556 << std::endl;
557
558 // Throw an error
562 }
563 } // if (dynamic_cast<GMRESBlockPreconditioner*>...
564 } // for (unsigned i=0;i<n_block;i++)
565
566 // Loop over the entries of the iteration count vector
567 for (unsigned i = 0; i < n_block; i++)
568 {
569 // If we used GMRES on this block
570 if (dynamic_cast<GMRESBlockPreconditioner*>(
571 this->Subsidiary_preconditioner_pt[i]) != 0)
572 {
573 // Update the mean
574 n_iter_mean +=
576 } // if (dynamic_cast<GMRESBlockPreconditioner*>...
577 } // for (unsigned i=0;i<n_block;i++)
578
579 // Loop over the entries of the iteration count vector (again)
580 // NOTE: We don't need to do the error checks again, we did them on the
581 // previous run through
582 for (unsigned i = 0; i < n_block; i++)
583 {
584 // If we used GMRES on this block
585 if (dynamic_cast<GMRESBlockPreconditioner*>(
586 this->Subsidiary_preconditioner_pt[i]) != 0)
587 {
588 // Update the variance
589 n_iter_var +=
590 (std::pow(double(iteration_count[i]) - n_iter_mean, 2.0) /
592 } // if (dynamic_cast<GMRESBlockPreconditioner*>...
593 } // for (unsigned i=0;i<n_block;i++)
594
595 // Doc. the statistics
596 oomph_info << "\nNumber of subsidiary blocks solved with GMRES: "
598 << "\nAverage subsidiary preconditioner iteration count: "
599 << n_iter_mean
600 << "\nStandard deviation of subsidiary preconditioner "
601 << "iteration count: " << std::sqrt(n_iter_var) << std::endl;
602 } // if (n_block_solved_with_gmres>0)
603 } // End of preconditioner_solve
604
605 // Ensure build of required objects (BUT only for CRDoubleMatrix objects)
608} // End of namespace oomph
cstr elem_len * i
Definition cfortran.h:603
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r.
A class for compressed row matrices. This is a distributable object.
Definition matrices.h:888
A vector in the mathematical sense, initially developed for linear algebra type applications....
void preconditioner_solve(const DoubleVector &r, DoubleVector &z)
Apply preconditioner to r.
The block preconditioner form of GMRES. This version extracts the blocks from the global systems and ...
Matrix vector product helper class - primarily a wrapper to Trilinos's Epetra matrix vector product m...
std::ostream *& stream_pt()
Access function for the stream pointer.
An OomphLibError object which should be thrown when an run-time error is encountered....
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
TAdvectionDiffusionReactionElement()
Constructor: Call constructors for TElement and AdvectionDiffusionReaction equations.
void split(const DoubleVector &in_vector, Vector< DoubleVector * > &out_vector_pt)
Split a DoubleVector into the out DoubleVectors. Let vec_A be the in Vector, and let vec_B and vec_C ...
void concatenate(const Vector< DoubleVector * > &in_vector_pt, DoubleVector &out_vector)
Concatenate DoubleVectors. Takes a Vector of DoubleVectors. If the out vector is built,...
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...