Reference documentation for deal.II version 8.5.1
petsc_precondition.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #include <deal.II/lac/petsc_precondition.h>
17 
18 #ifdef DEAL_II_WITH_PETSC
19 
20 # include <deal.II/base/utilities.h>
21 # include <deal.II/lac/exceptions.h>
22 # include <deal.II/lac/petsc_compatibility.h>
23 # include <deal.II/lac/petsc_matrix_base.h>
24 # include <deal.II/lac/petsc_vector_base.h>
25 # include <deal.II/lac/petsc_solver.h>
26 # include <petscconf.h>
27 # include <cmath>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 namespace PETScWrappers
32 {
34  :
35  pc(NULL), matrix(NULL)
36  {}
37 
39  {
40  try
41  {
42  clear();
43  }
44  catch (...)
45  {}
46  }
47 
48  void
50  {
51  matrix = NULL;
52 
53  if (pc!=NULL)
54  {
55 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
56  PetscErrorCode ierr = PCDestroy(pc);
57 #else
58  PetscErrorCode ierr = PCDestroy(&pc);
59 #endif
60  pc = NULL;
61  AssertThrow (ierr == 0, ExcPETScError(ierr));
62  }
63  }
64 
65 
66  void
68  const VectorBase &src) const
69  {
71 
72  const PetscErrorCode ierr = PCApply(pc, src, dst);
73  AssertThrow (ierr == 0, ExcPETScError(ierr));
74  }
75 
76 
77  void
79  {
80  // only allow the creation of the
81  // preconditioner once
83 
84  MPI_Comm comm;
85  // this ugly cast is necessary because the
86  // type Mat and PETScObject are
87  // unrelated.
88  PetscErrorCode ierr = PetscObjectGetComm(reinterpret_cast<PetscObject>(matrix),
89  &comm);
90  AssertThrow (ierr == 0, ExcPETScError(ierr));
91 
92  ierr = PCCreate(comm, &pc);
93  AssertThrow (ierr == 0, ExcPETScError(ierr));
94 
95 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
96  ierr = PCSetOperators(pc , matrix, matrix, SAME_PRECONDITIONER);
97 #else
98  ierr = PCSetOperators(pc , matrix, matrix);
99 #endif
100  AssertThrow (ierr == 0, ExcPETScError(ierr));
101  }
102 
103 
104  const PC &
106  {
107  return pc;
108  }
109 
110 
111  PreconditionerBase::operator Mat () const
112  {
113  return matrix;
114  }
115 
116 
117  /* ----------------- PreconditionJacobi -------------------- */
119  const AdditionalData &additional_data_)
120  {
121  additional_data = additional_data_;
122 
123  PetscErrorCode ierr = PCCreate(comm, &pc);
124  AssertThrow (ierr == 0, ExcPETScError(ierr));
125 
126  initialize();
127  }
128 
129 
131  {}
132 
133 
135  const AdditionalData &additional_data)
136  {
138  }
139 
140  void
142  {
144 
145  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCJACOBI));
146  AssertThrow (ierr == 0, ExcPETScError(ierr));
147 
148  ierr = PCSetFromOptions (pc);
149  AssertThrow (ierr == 0, ExcPETScError(ierr));
150  }
151 
152  void
154  const AdditionalData &additional_data_)
155  {
156  clear ();
157 
158  matrix = static_cast<Mat>(matrix_);
159  additional_data = additional_data_;
160 
161  create_pc();
162  initialize();
163 
164  PetscErrorCode ierr = PCSetUp (pc);
165  AssertThrow (ierr == 0, ExcPETScError(ierr));
166  }
167 
168 
169  /* ----------------- PreconditionBlockJacobi -------------------- */
171  const AdditionalData &additional_data_)
172  {
173  additional_data = additional_data_;
174 
175  PetscErrorCode ierr = PCCreate(comm, &pc);
176  AssertThrow (ierr == 0, ExcPETScError(ierr));
177 
178  initialize();
179  }
180 
181 
183  {}
184 
185 
188  const AdditionalData &additional_data)
189  {
191  }
192 
193  void
195  {
196  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCBJACOBI));
197  AssertThrow (ierr == 0, ExcPETScError(ierr));
198 
199  ierr = PCSetFromOptions (pc);
200  AssertThrow (ierr == 0, ExcPETScError(ierr));
201  }
202 
203 
204  void
206  const AdditionalData &additional_data_)
207  {
208  clear ();
209 
210  matrix = static_cast<Mat>(matrix_);
211  additional_data = additional_data_;
212 
213  create_pc();
214  initialize();
215 
216  PetscErrorCode ierr = PCSetUp (pc);
217  AssertThrow (ierr == 0, ExcPETScError(ierr));
218  }
219 
220 
221  /* ----------------- PreconditionSOR -------------------- */
222 
224  AdditionalData (const double omega)
225  :
226  omega (omega)
227  {}
228 
229 
231  {}
232 
233 
236  {
238  }
239 
240 
241  void
243  const AdditionalData &additional_data_)
244  {
245  clear ();
246 
247  matrix = static_cast<Mat>(matrix_);
248  additional_data = additional_data_;
249 
250  create_pc();
251 
252  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCSOR));
253  AssertThrow (ierr == 0, ExcPETScError(ierr));
254 
255  // then set flags as given
256  ierr = PCSORSetOmega (pc, additional_data.omega);
257  AssertThrow (ierr == 0, ExcPETScError(ierr));
258 
259  ierr = PCSetFromOptions (pc);
260  AssertThrow (ierr == 0, ExcPETScError(ierr));
261 
262  ierr = PCSetUp (pc);
263  AssertThrow (ierr == 0, ExcPETScError(ierr));
264  }
265 
266 
267  /* ----------------- PreconditionSSOR -------------------- */
268 
270  AdditionalData (const double omega)
271  :
272  omega (omega)
273  {}
274 
275 
277  {}
278 
279 
282  {
284  }
285 
286 
287  void
289  const AdditionalData &additional_data_)
290  {
291  clear ();
292 
293  matrix = static_cast<Mat>(matrix_);
294  additional_data = additional_data_;
295 
296  create_pc();
297 
298  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCSOR));
299  AssertThrow (ierr == 0, ExcPETScError(ierr));
300 
301  // then set flags as given
302  ierr = PCSORSetOmega (pc, additional_data.omega);
303  AssertThrow (ierr == 0, ExcPETScError(ierr));
304 
305  // convert SOR to SSOR
306  ierr = PCSORSetSymmetric (pc, SOR_SYMMETRIC_SWEEP);
307  AssertThrow (ierr == 0, ExcPETScError(ierr));
308 
309  ierr = PCSetFromOptions (pc);
310  AssertThrow (ierr == 0, ExcPETScError(ierr));
311 
312  ierr = PCSetUp (pc);
313  AssertThrow (ierr == 0, ExcPETScError(ierr));
314  }
315 
316 
317  /* ----------------- PreconditionEisenstat -------------------- */
318 
320  AdditionalData (const double omega)
321  :
322  omega (omega)
323  {}
324 
325 
327  {}
328 
329 
332  {
334  }
335 
336 
337  void
339  const AdditionalData &additional_data_)
340  {
341  clear ();
342 
343  matrix = static_cast<Mat>(matrix_);
344  additional_data = additional_data_;
345 
346  create_pc();
347 
348  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCEISENSTAT));
349  AssertThrow (ierr == 0, ExcPETScError(ierr));
350 
351  // then set flags as given
352  ierr = PCEisenstatSetOmega (pc, additional_data.omega);
353  AssertThrow (ierr == 0, ExcPETScError(ierr));
354 
355  ierr = PCSetFromOptions (pc);
356  AssertThrow (ierr == 0, ExcPETScError(ierr));
357 
358  ierr = PCSetUp (pc);
359  AssertThrow (ierr == 0, ExcPETScError(ierr));
360  }
361 
362 
363  /* ----------------- PreconditionICC -------------------- */
364 
365 
367  AdditionalData (const unsigned int levels)
368  :
369  levels (levels)
370  {}
371 
372 
374  {}
375 
376 
379  {
381  }
382 
383 
384  void
386  const AdditionalData &additional_data_)
387  {
388  clear ();
389 
390  matrix = static_cast<Mat>(matrix_);
391  additional_data = additional_data_;
392 
393  create_pc();
394 
395  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCICC));
396  AssertThrow (ierr == 0, ExcPETScError(ierr));
397 
398  // then set flags
399  ierr = PCFactorSetLevels (pc, additional_data.levels);
400  AssertThrow (ierr == 0, ExcPETScError(ierr));
401 
402  ierr = PCSetFromOptions (pc);
403  AssertThrow (ierr == 0, ExcPETScError(ierr));
404 
405  ierr = PCSetUp (pc);
406  AssertThrow (ierr == 0, ExcPETScError(ierr));
407  }
408 
409 
410  /* ----------------- PreconditionILU -------------------- */
411 
413  AdditionalData (const unsigned int levels)
414  :
415  levels (levels)
416  {}
417 
418 
420  {}
421 
422 
425  {
427  }
428 
429 
430  void
432  const AdditionalData &additional_data_)
433  {
434  clear ();
435 
436  matrix = static_cast<Mat>(matrix_);
437  additional_data = additional_data_;
438 
439  create_pc();
440 
441  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCILU));
442  AssertThrow (ierr == 0, ExcPETScError(ierr));
443 
444  // then set flags
445  ierr = PCFactorSetLevels (pc, additional_data.levels);
446  AssertThrow (ierr == 0, ExcPETScError(ierr));
447 
448  ierr = PCSetFromOptions (pc);
449  AssertThrow (ierr == 0, ExcPETScError(ierr));
450 
451  ierr = PCSetUp (pc);
452  AssertThrow (ierr == 0, ExcPETScError(ierr));
453  }
454 
455 
456  /* ----------------- PreconditionBoomerAMG -------------------- */
457 
459  AdditionalData(const bool symmetric_operator,
460  const double strong_threshold,
461  const double max_row_sum,
462  const unsigned int aggressive_coarsening_num_levels,
463  const bool output_details
464  )
465  :
466  symmetric_operator(symmetric_operator),
467  strong_threshold(strong_threshold),
468  max_row_sum(max_row_sum),
469  aggressive_coarsening_num_levels(aggressive_coarsening_num_levels),
470  output_details(output_details)
471  {}
472 
473 
475  {}
476 
478  const AdditionalData &additional_data_)
479  {
480  additional_data = additional_data_;
481 
482  PetscErrorCode ierr = PCCreate(comm, &pc);
483  AssertThrow (ierr == 0, ExcPETScError(ierr));
484 
485 #ifdef PETSC_HAVE_HYPRE
486  initialize();
487 #else // PETSC_HAVE_HYPRE
488  (void)pc;
489  Assert (false,
490  ExcMessage ("Your PETSc installation does not include a copy of "
491  "the hypre package necessary for this preconditioner."));
492 #endif
493  }
494 
495 
498  {
500  }
501 
502  void
504  {
505 #ifdef PETSC_HAVE_HYPRE
506  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCHYPRE));
507  AssertThrow (ierr == 0, ExcPETScError(ierr));
508 
509  ierr = PCHYPRESetType(pc, "boomeramg");
510  AssertThrow (ierr == 0, ExcPETScError(ierr));
511 
513  {
514  set_option_value("-pc_hypre_boomeramg_print_statistics", "1");
515  }
516 
517  set_option_value("-pc_hypre_boomeramg_agg_nl",
519 
520  std::stringstream ssStream;
521  ssStream << additional_data.max_row_sum;
522  set_option_value("-pc_hypre_boomeramg_max_row_sum", ssStream.str());
523 
524  ssStream.str(""); // empty the stringstream
525  ssStream << additional_data.strong_threshold;
526  set_option_value("-pc_hypre_boomeramg_strong_threshold", ssStream.str());
527 
529  {
530  set_option_value("-pc_hypre_boomeramg_relax_type_up", "symmetric-SOR/Jacobi");
531  set_option_value("-pc_hypre_boomeramg_relax_type_down", "symmetric-SOR/Jacobi");
532  set_option_value("-pc_hypre_boomeramg_relax_type_coarse", "Gaussian-elimination");
533  }
534  else
535  {
536  set_option_value("-pc_hypre_boomeramg_relax_type_up", "SOR/Jacobi");
537  set_option_value("-pc_hypre_boomeramg_relax_type_down", "SOR/Jacobi");
538  set_option_value("-pc_hypre_boomeramg_relax_type_coarse", "Gaussian-elimination");
539  }
540 
541  ierr = PCSetFromOptions (pc);
542  AssertThrow (ierr == 0, ExcPETScError(ierr));
543 #else
544  Assert (false,
545  ExcMessage ("Your PETSc installation does not include a copy of "
546  "the hypre package necessary for this preconditioner."));
547 #endif
548  }
549 
550  void
552  const AdditionalData &additional_data_)
553  {
554 #ifdef PETSC_HAVE_HYPRE
555  clear ();
556 
557  matrix = static_cast<Mat>(matrix_);
558  additional_data = additional_data_;
559 
560  create_pc();
561  initialize ();
562 
563  PetscErrorCode ierr = PCSetUp (pc);
564  AssertThrow (ierr == 0, ExcPETScError(ierr));
565 
566 #else // PETSC_HAVE_HYPRE
567  (void)matrix_;
568  (void)additional_data_;
569  Assert (false,
570  ExcMessage ("Your PETSc installation does not include a copy of "
571  "the hypre package necessary for this preconditioner."));
572 #endif
573  }
574 
575 
576  /* ----------------- PreconditionParaSails -------------------- */
577 
579  AdditionalData(const unsigned int symmetric,
580  const unsigned int n_levels,
581  const double threshold,
582  const double filter,
583  const bool output_details)
584  :
585  symmetric(symmetric),
586  n_levels(n_levels),
587  threshold(threshold),
588  filter(filter),
589  output_details(output_details)
590  {}
591 
592 
594  {}
595 
596 
599  {
601  }
602 
603 
604  void
606  const AdditionalData &additional_data_)
607  {
608  clear ();
609 
610  matrix = static_cast<Mat>(matrix_);
611  additional_data = additional_data_;
612 
613 #ifdef PETSC_HAVE_HYPRE
614  create_pc();
615 
616  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCHYPRE));
617  AssertThrow (ierr == 0, ExcPETScError(ierr));
618 
619  ierr = PCHYPRESetType(pc, "parasails");
620  AssertThrow (ierr == 0, ExcPETScError(ierr));
621 
623  {
624  set_option_value("-pc_hypre_parasails_logging","1");
625  }
626 
628  additional_data.symmetric == 1 ||
630  ExcMessage("ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
631 
632  std::stringstream ssStream;
633 
634  switch (additional_data.symmetric)
635  {
636  case 0:
637  {
638  ssStream << "nonsymmetric";
639  break;
640  }
641 
642  case 1:
643  {
644  ssStream << "SPD";
645  break;
646  }
647 
648  case 2:
649  {
650  ssStream << "nonsymmetric,SPD";
651  break;
652  }
653 
654  default:
655  Assert (false,
656  ExcMessage("ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
657  };
658 
659  set_option_value("-pc_hypre_parasails_sym",ssStream.str());
660 
661  set_option_value ("-pc_hypre_parasails_nlevels",
663 
664  ssStream.str(""); // empty the stringstream
665  ssStream << additional_data.threshold;
666  set_option_value("-pc_hypre_parasails_thresh", ssStream.str());
667 
668  ssStream.str(""); // empty the stringstream
669  ssStream << additional_data.filter;
670  set_option_value("-pc_hypre_parasails_filter", ssStream.str());
671 
672  ierr = PCSetFromOptions (pc);
673  AssertThrow (ierr == 0, ExcPETScError(ierr));
674 
675  ierr = PCSetUp (pc);
676  AssertThrow (ierr == 0, ExcPETScError(ierr));
677 
678 #else // PETSC_HAVE_HYPRE
679  (void)pc;
680  Assert (false,
681  ExcMessage ("Your PETSc installation does not include a copy of "
682  "the hypre package necessary for this preconditioner."));
683 #endif
684  }
685 
686 
687  /* ----------------- PreconditionNone ------------------------- */
688 
690  {}
691 
692 
695  {
697  }
698 
699 
700  void
702  const AdditionalData &additional_data_)
703  {
704  clear ();
705 
706  matrix = static_cast<Mat>(matrix_);
707  additional_data = additional_data_;
708 
709  create_pc();
710 
711  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCNONE));
712  AssertThrow (ierr == 0, ExcPETScError(ierr));
713 
714  ierr = PCSetFromOptions (pc);
715  AssertThrow (ierr == 0, ExcPETScError(ierr));
716 
717  ierr = PCSetUp (pc);
718  AssertThrow (ierr == 0, ExcPETScError(ierr));
719  }
720 
721 
722  /* ----------------- PreconditionLU -------------------- */
723 
725  AdditionalData (const double pivoting,
726  const double zero_pivot,
727  const double damping)
728  :
729  pivoting (pivoting),
730  zero_pivot (zero_pivot),
731  damping (damping)
732  {}
733 
734 
736  {}
737 
738 
741  {
743  }
744 
745 
746  void
748  const AdditionalData &additional_data_)
749  {
750  clear ();
751 
752  matrix = static_cast<Mat>(matrix_);
753  additional_data = additional_data_;
754 
755  create_pc();
756 
757  PetscErrorCode ierr = PCSetType (pc, const_cast<char *>(PCLU));
758  AssertThrow (ierr == 0, ExcPETScError(ierr));
759 
760  // set flags as given
761 #if DEAL_II_PETSC_VERSION_LT(3,0,1)
762  ierr = PCFactorSetPivoting (pc, additional_data.pivoting);
763 #else
764  ierr = PCFactorSetColumnPivot (pc, additional_data.pivoting);
765 #endif
766  AssertThrow (ierr == 0, ExcPETScError(ierr));
767 
768  ierr = PCFactorSetZeroPivot (pc, additional_data.zero_pivot);
769  AssertThrow (ierr == 0, ExcPETScError(ierr));
770 
771 #if DEAL_II_PETSC_VERSION_LT(3,0,1)
772  ierr = PCFactorSetShiftNonzero (pc, additional_data.damping);
773 #else
774  ierr = PCFactorSetShiftAmount (pc, additional_data.damping);
775 #endif
776  AssertThrow (ierr == 0, ExcPETScError(ierr));
777 
778  ierr = PCSetFromOptions (pc);
779  AssertThrow (ierr == 0, ExcPETScError(ierr));
780 
781  ierr = PCSetUp (pc);
782  AssertThrow (ierr == 0, ExcPETScError(ierr));
783  }
784 
785 }
786 
787 DEAL_II_NAMESPACE_CLOSE
788 
789 #endif // DEAL_II_WITH_PETSC
void vmult(VectorBase &dst, const VectorBase &src) const
AdditionalData(const double pivoting=1.e-6, const double zero_pivot=1.e-12, const double damping=0.0)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
#define AssertThrow(cond, exc)
Definition: exceptions.h:369
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:92
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void set_option_value(const std::string &name, const std::string &value)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
static ::ExceptionBase & ExcInvalidState()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const unsigned int symmetric=1, const unsigned int n_levels=1, const double threshold=0.1, const double filter=0.05, const bool output_details=false)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
AdditionalData(const bool symmetric_operator=false, const double strong_threshold=0.25, const double max_row_sum=0.9, const unsigned int aggressive_coarsening_num_levels=0, const bool output_details=false)