Reference documentation for deal.II version 8.4.2
petsc_precondition.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2015 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/petsc_compatibility.h>
22 # include <deal.II/lac/petsc_matrix_base.h>
23 # include <deal.II/lac/petsc_vector_base.h>
24 # include <deal.II/lac/petsc_solver.h>
25 # include <petscconf.h>
26 # include <cmath>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 namespace PETScWrappers
31 {
33  :
34  pc(NULL), matrix(NULL)
35  {}
36 
37 
39  {
40  if (pc!=NULL)
41  {
42 #if DEAL_II_PETSC_VERSION_LT(3,2,0)
43  int ierr = PCDestroy(pc);
44 #else
45  int ierr = PCDestroy(&pc);
46 #endif
47  AssertThrow (ierr == 0, ExcPETScError(ierr));
48  }
49  }
50 
51 
52  void
54  const VectorBase &src) const
55  {
56  AssertThrow (pc != NULL, StandardExceptions::ExcInvalidState ());
57 
58  int ierr;
59  ierr = PCApply(pc, src, dst);
60  AssertThrow (ierr == 0, ExcPETScError(ierr));
61  }
62 
63 
64  void
66  {
67  // only allow the creation of the
68  // preconditioner once
69  AssertThrow (pc == NULL, StandardExceptions::ExcInvalidState ());
70 
71  MPI_Comm comm;
72  int ierr;
73  // this ugly cast is necessary because the
74  // type Mat and PETScObject are
75  // unrelated.
76  ierr = PetscObjectGetComm(reinterpret_cast<PetscObject>(matrix), &comm);
77  AssertThrow (ierr == 0, ExcPETScError(ierr));
78 
79  ierr = PCCreate(comm, &pc);
80  AssertThrow (ierr == 0, ExcPETScError(ierr));
81 
82 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
83  ierr = PCSetOperators(pc , matrix, matrix, SAME_PRECONDITIONER);
84 #else
85  ierr = PCSetOperators(pc , matrix, matrix);
86 #endif
87  AssertThrow (ierr == 0, ExcPETScError(ierr));
88  }
89 
90 
91  const PC &
93  {
94  return pc;
95  }
96 
97 
98  PreconditionerBase::operator Mat () const
99  {
100  return matrix;
101  }
102 
103 
104  /* ----------------- PreconditionJacobi -------------------- */
106  const AdditionalData &additional_data_)
107  {
108  additional_data = additional_data_;
109 
110  int ierr = PCCreate(comm, &pc);
111  AssertThrow (ierr == 0, ExcPETScError(ierr));
112 
113  initialize();
114  }
115 
116 
118  {}
119 
120 
122  const AdditionalData &additional_data)
123  {
124  initialize(matrix, additional_data);
125  }
126 
127  void
129  {
130  int ierr;
131  ierr = PCSetType (pc, const_cast<char *>(PCJACOBI));
132  AssertThrow (ierr == 0, ExcPETScError(ierr));
133 
134  ierr = PCSetFromOptions (pc);
135  AssertThrow (ierr == 0, ExcPETScError(ierr));
136  }
137 
138  void
140  const AdditionalData &additional_data_)
141  {
142  matrix = static_cast<Mat>(matrix_);
143  additional_data = additional_data_;
144 
145  create_pc();
146  initialize();
147 
148  int ierr = PCSetUp (pc);
149  AssertThrow (ierr == 0, ExcPETScError(ierr));
150  }
151 
152 
153  /* ----------------- PreconditionBlockJacobi -------------------- */
155  const AdditionalData &additional_data_)
156  {
157  additional_data = additional_data_;
158 
159  int ierr = PCCreate(comm, &pc);
160  AssertThrow (ierr == 0, ExcPETScError(ierr));
161 
162  initialize();
163  }
164 
165 
167  {}
168 
169 
172  const AdditionalData &additional_data)
173  {
174  initialize(matrix, additional_data);
175  }
176 
177  void
179  {
180  int ierr;
181  ierr = PCSetType (pc, const_cast<char *>(PCBJACOBI));
182  AssertThrow (ierr == 0, ExcPETScError(ierr));
183 
184  ierr = PCSetFromOptions (pc);
185  AssertThrow (ierr == 0, ExcPETScError(ierr));
186  }
187 
188 
189  void
191  const AdditionalData &additional_data_)
192  {
193  matrix = static_cast<Mat>(matrix_);
194  additional_data = additional_data_;
195 
196  create_pc();
197  initialize();
198 
199  int ierr = PCSetUp (pc);
200  AssertThrow (ierr == 0, ExcPETScError(ierr));
201  }
202 
203 
204  /* ----------------- PreconditionSOR -------------------- */
205 
207  AdditionalData (const double omega)
208  :
209  omega (omega)
210  {}
211 
212 
214  {}
215 
216 
219  {
220  initialize(matrix, additional_data);
221  }
222 
223 
224  void
226  const AdditionalData &additional_data_)
227  {
228  matrix = static_cast<Mat>(matrix_);
229  additional_data = additional_data_;
230 
231  create_pc();
232 
233  int ierr;
234  ierr = PCSetType (pc, const_cast<char *>(PCSOR));
235  AssertThrow (ierr == 0, ExcPETScError(ierr));
236 
237  // then set flags as given
238  ierr = PCSORSetOmega (pc, additional_data.omega);
239  AssertThrow (ierr == 0, ExcPETScError(ierr));
240 
241  ierr = PCSetFromOptions (pc);
242  AssertThrow (ierr == 0, ExcPETScError(ierr));
243 
244  ierr = PCSetUp (pc);
245  AssertThrow (ierr == 0, ExcPETScError(ierr));
246  }
247 
248 
249  /* ----------------- PreconditionSSOR -------------------- */
250 
252  AdditionalData (const double omega)
253  :
254  omega (omega)
255  {}
256 
257 
259  {}
260 
261 
264  {
265  initialize(matrix, additional_data);
266  }
267 
268 
269  void
271  const AdditionalData &additional_data_)
272  {
273  matrix = static_cast<Mat>(matrix_);
274  additional_data = additional_data_;
275 
276  create_pc();
277 
278  int ierr;
279  ierr = PCSetType (pc, const_cast<char *>(PCSOR));
280  AssertThrow (ierr == 0, ExcPETScError(ierr));
281 
282  // then set flags as given
283  ierr = PCSORSetOmega (pc, additional_data.omega);
284  AssertThrow (ierr == 0, ExcPETScError(ierr));
285 
286  // convert SOR to SSOR
287  ierr = PCSORSetSymmetric (pc, SOR_SYMMETRIC_SWEEP);
288  AssertThrow (ierr == 0, ExcPETScError(ierr));
289 
290  ierr = PCSetFromOptions (pc);
291  AssertThrow (ierr == 0, ExcPETScError(ierr));
292 
293  ierr = PCSetUp (pc);
294  AssertThrow (ierr == 0, ExcPETScError(ierr));
295  }
296 
297 
298  /* ----------------- PreconditionEisenstat -------------------- */
299 
301  AdditionalData (const double omega)
302  :
303  omega (omega)
304  {}
305 
306 
308  {}
309 
310 
313  {
314  initialize(matrix, additional_data);
315  }
316 
317 
318  void
320  const AdditionalData &additional_data_)
321  {
322  matrix = static_cast<Mat>(matrix_);
323  additional_data = additional_data_;
324 
325  create_pc();
326 
327  int ierr;
328  ierr = PCSetType (pc, const_cast<char *>(PCEISENSTAT));
329  AssertThrow (ierr == 0, ExcPETScError(ierr));
330 
331  // then set flags as given
332  ierr = PCEisenstatSetOmega (pc, additional_data.omega);
333  AssertThrow (ierr == 0, ExcPETScError(ierr));
334 
335  ierr = PCSetFromOptions (pc);
336  AssertThrow (ierr == 0, ExcPETScError(ierr));
337 
338  ierr = PCSetUp (pc);
339  AssertThrow (ierr == 0, ExcPETScError(ierr));
340  }
341 
342 
343  /* ----------------- PreconditionICC -------------------- */
344 
345 
347  AdditionalData (const unsigned int levels)
348  :
349  levels (levels)
350  {}
351 
352 
354  {}
355 
356 
359  {
360  initialize(matrix, additional_data);
361  }
362 
363 
364  void
366  const AdditionalData &additional_data_)
367  {
368  matrix = static_cast<Mat>(matrix_);
369  additional_data = additional_data_;
370 
371  create_pc();
372 
373  int ierr;
374  ierr = PCSetType (pc, const_cast<char *>(PCICC));
375  AssertThrow (ierr == 0, ExcPETScError(ierr));
376 
377  // then set flags
378  PCFactorSetLevels (pc, additional_data.levels);
379  AssertThrow (ierr == 0, ExcPETScError(ierr));
380 
381  ierr = PCSetFromOptions (pc);
382  AssertThrow (ierr == 0, ExcPETScError(ierr));
383 
384  ierr = PCSetUp (pc);
385  AssertThrow (ierr == 0, ExcPETScError(ierr));
386  }
387 
388 
389  /* ----------------- PreconditionILU -------------------- */
390 
392  AdditionalData (const unsigned int levels)
393  :
394  levels (levels)
395  {}
396 
397 
399  {}
400 
401 
404  {
405  initialize(matrix, additional_data);
406  }
407 
408 
409  void
411  const AdditionalData &additional_data_)
412  {
413  matrix = static_cast<Mat>(matrix_);
414  additional_data = additional_data_;
415 
416  create_pc();
417 
418  int ierr;
419  ierr = PCSetType (pc, const_cast<char *>(PCILU));
420  AssertThrow (ierr == 0, ExcPETScError(ierr));
421 
422  // then set flags
423  PCFactorSetLevels (pc, additional_data.levels);
424  AssertThrow (ierr == 0, ExcPETScError(ierr));
425 
426  ierr = PCSetFromOptions (pc);
427  AssertThrow (ierr == 0, ExcPETScError(ierr));
428 
429  ierr = PCSetUp (pc);
430  AssertThrow (ierr == 0, ExcPETScError(ierr));
431  }
432 
433 
434  /* ----------------- PreconditionBoomerAMG -------------------- */
435 
437  AdditionalData(const bool symmetric_operator,
438  const double strong_threshold,
439  const double max_row_sum,
440  const unsigned int aggressive_coarsening_num_levels,
441  const bool output_details
442  )
443  :
444  symmetric_operator(symmetric_operator),
445  strong_threshold(strong_threshold),
446  max_row_sum(max_row_sum),
447  aggressive_coarsening_num_levels(aggressive_coarsening_num_levels),
448  output_details(output_details)
449  {}
450 
451 
453  {}
454 
456  const AdditionalData &additional_data_)
457  {
458  additional_data = additional_data_;
459 
460  int ierr = PCCreate(comm, &pc);
461  AssertThrow (ierr == 0, ExcPETScError(ierr));
462 
463 #ifdef PETSC_HAVE_HYPRE
464  initialize();
465 #else // PETSC_HAVE_HYPRE
466  (void)pc;
467  Assert (false,
468  ExcMessage ("Your PETSc installation does not include a copy of "
469  "the hypre package necessary for this preconditioner."));
470 #endif
471  }
472 
473 
476  {
477  initialize(matrix, additional_data);
478  }
479 
480  void
482  {
483  int ierr;
484  ierr = PCSetType (pc, const_cast<char *>(PCHYPRE));
485  AssertThrow (ierr == 0, ExcPETScError(ierr));
486 
487  ierr = PCHYPRESetType(pc, "boomeramg");
488  AssertThrow (ierr == 0, ExcPETScError(ierr));
489 
491  {
492  set_option_value("-pc_hypre_boomeramg_print_statistics", "1");
493  }
494 
495  set_option_value("-pc_hypre_boomeramg_agg_nl",
497 
498  std::stringstream ssStream;
499  ssStream << additional_data.max_row_sum;
500  set_option_value("-pc_hypre_boomeramg_max_row_sum", ssStream.str());
501 
502  ssStream.str(""); // empty the stringstream
503  ssStream << additional_data.strong_threshold;
504  set_option_value("-pc_hypre_boomeramg_strong_threshold", ssStream.str());
505 
507  {
508  set_option_value("-pc_hypre_boomeramg_relax_type_up", "symmetric-SOR/Jacobi");
509  set_option_value("-pc_hypre_boomeramg_relax_type_down", "symmetric-SOR/Jacobi");
510  set_option_value("-pc_hypre_boomeramg_relax_type_coarse", "Gaussian-elimination");
511  }
512  else
513  {
514  set_option_value("-pc_hypre_boomeramg_relax_type_up", "SOR/Jacobi");
515  set_option_value("-pc_hypre_boomeramg_relax_type_down", "SOR/Jacobi");
516  set_option_value("-pc_hypre_boomeramg_relax_type_coarse", "Gaussian-elimination");
517  }
518 
519  ierr = PCSetFromOptions (pc);
520  AssertThrow (ierr == 0, ExcPETScError(ierr));
521  }
522 
523  void
525  const AdditionalData &additional_data_)
526  {
527  matrix = static_cast<Mat>(matrix_);
528  additional_data = additional_data_;
529 
530 #ifdef PETSC_HAVE_HYPRE
531  create_pc();
532  initialize ();
533 
534  int ierr = PCSetUp (pc);
535  AssertThrow (ierr == 0, ExcPETScError(ierr));
536 
537 #else // PETSC_HAVE_HYPRE
538  (void)pc;
539  Assert (false,
540  ExcMessage ("Your PETSc installation does not include a copy of "
541  "the hypre package necessary for this preconditioner."));
542 #endif
543  }
544 
545 
546  /* ----------------- PreconditionParaSails -------------------- */
547 
549  AdditionalData(const unsigned int symmetric,
550  const unsigned int n_levels,
551  const double threshold,
552  const double filter,
553  const bool output_details)
554  :
555  symmetric(symmetric),
556  n_levels(n_levels),
557  threshold(threshold),
558  filter(filter),
559  output_details(output_details)
560  {}
561 
562 
564  {}
565 
566 
569  {
570  initialize(matrix, additional_data);
571  }
572 
573 
574  void
576  const AdditionalData &additional_data_)
577  {
578  matrix = static_cast<Mat>(matrix_);
579  additional_data = additional_data_;
580 
581 #ifdef PETSC_HAVE_HYPRE
582  create_pc();
583 
584  int ierr;
585  ierr = PCSetType (pc, const_cast<char *>(PCHYPRE));
586  AssertThrow (ierr == 0, ExcPETScError(ierr));
587 
588  ierr = PCHYPRESetType(pc, "parasails");
589  AssertThrow (ierr == 0, ExcPETScError(ierr));
590 
592  {
593  set_option_value("-pc_hypre_parasails_logging","1");
594  }
595 
597  additional_data.symmetric == 1 ||
599  ExcMessage("ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
600 
601  std::stringstream ssStream;
602 
603  switch (additional_data.symmetric)
604  {
605  case 0:
606  {
607  ssStream << "nonsymmetric";
608  break;
609  }
610 
611  case 1:
612  {
613  ssStream << "SPD";
614  break;
615  }
616 
617  case 2:
618  {
619  ssStream << "nonsymmetric,SPD";
620  break;
621  }
622 
623  default:
624  Assert (false,
625  ExcMessage("ParaSails parameter symmetric can only be equal to 0, 1, 2!"));
626  };
627 
628  set_option_value("-pc_hypre_parasails_sym",ssStream.str());
629 
630  set_option_value ("-pc_hypre_parasails_nlevels",
632 
633  ssStream.str(""); // empty the stringstream
634  ssStream << additional_data.threshold;
635  set_option_value("-pc_hypre_parasails_thresh", ssStream.str());
636 
637  ssStream.str(""); // empty the stringstream
638  ssStream << additional_data.filter;
639  set_option_value("-pc_hypre_parasails_filter", ssStream.str());
640 
641  ierr = PCSetFromOptions (pc);
642  AssertThrow (ierr == 0, ExcPETScError(ierr));
643 
644  ierr = PCSetUp (pc);
645  AssertThrow (ierr == 0, ExcPETScError(ierr));
646 
647 #else // PETSC_HAVE_HYPRE
648  (void)pc;
649  Assert (false,
650  ExcMessage ("Your PETSc installation does not include a copy of "
651  "the hypre package necessary for this preconditioner."));
652 #endif
653  }
654 
655 
656  /* ----------------- PreconditionNone ------------------------- */
657 
659  {}
660 
661 
664  {
665  initialize (matrix, additional_data);
666  }
667 
668 
669  void
671  const AdditionalData &additional_data_)
672  {
673  matrix = static_cast<Mat>(matrix_);
674  additional_data = additional_data_;
675 
676  create_pc();
677 
678  int ierr;
679  ierr = PCSetType (pc, const_cast<char *>(PCNONE));
680  AssertThrow (ierr == 0, ExcPETScError(ierr));
681 
682  ierr = PCSetFromOptions (pc);
683  AssertThrow (ierr == 0, ExcPETScError(ierr));
684 
685  ierr = PCSetUp (pc);
686  AssertThrow (ierr == 0, ExcPETScError(ierr));
687  }
688 
689 
690  /* ----------------- PreconditionLU -------------------- */
691 
693  AdditionalData (const double pivoting,
694  const double zero_pivot,
695  const double damping)
696  :
697  pivoting (pivoting),
698  zero_pivot (zero_pivot),
699  damping (damping)
700  {}
701 
702 
704  {}
705 
706 
709  {
710  initialize(matrix, additional_data);
711  }
712 
713 
714  void
716  const AdditionalData &additional_data_)
717  {
718  matrix = static_cast<Mat>(matrix_);
719  additional_data = additional_data_;
720 
721  create_pc();
722 
723  int ierr;
724  ierr = PCSetType (pc, const_cast<char *>(PCLU));
725  AssertThrow (ierr == 0, ExcPETScError(ierr));
726 
727  // set flags as given
728 #if DEAL_II_PETSC_VERSION_LT(3,0,1)
729  ierr = PCFactorSetPivoting (pc, additional_data.pivoting);
730 #else
731  ierr = PCFactorSetColumnPivot (pc, additional_data.pivoting);
732 #endif
733  AssertThrow (ierr == 0, ExcPETScError(ierr));
734 
735  ierr = PCFactorSetZeroPivot (pc, additional_data.zero_pivot);
736  AssertThrow (ierr == 0, ExcPETScError(ierr));
737 
738 #if DEAL_II_PETSC_VERSION_LT(3,0,1)
739  ierr = PCFactorSetShiftNonzero (pc, additional_data.damping);
740 #else
741  ierr = PCFactorSetShiftAmount (pc, additional_data.damping);
742 #endif
743  AssertThrow (ierr == 0, ExcPETScError(ierr));
744 
745  ierr = PCSetFromOptions (pc);
746  AssertThrow (ierr == 0, ExcPETScError(ierr));
747 
748  ierr = PCSetUp (pc);
749  AssertThrow (ierr == 0, ExcPETScError(ierr));
750  }
751 
752 }
753 
754 DEAL_II_NAMESPACE_CLOSE
755 
756 #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())
::ExceptionBase & ExcMessage(std::string arg1)
void initialize(const MatrixBase &matrix, const AdditionalData &additional_data=AdditionalData())
#define AssertThrow(cond, exc)
Definition: exceptions.h:358
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:89
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())
#define Assert(cond, exc)
Definition: exceptions.h:294
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)