Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
petsc_solver.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2004 - 2018 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.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 
17 #include <deal.II/base/logstream.h>
18 #include <deal.II/base/std_cxx14/memory.h>
19 
20 #include <deal.II/lac/petsc_solver.h>
21 
22 #ifdef DEAL_II_WITH_PETSC
23 
24 # include <deal.II/lac/exceptions.h>
25 # include <deal.II/lac/petsc_compatibility.h>
26 # include <deal.II/lac/petsc_matrix_base.h>
27 # include <deal.II/lac/petsc_precondition.h>
28 # include <deal.II/lac/petsc_vector_base.h>
29 
30 # include <petscversion.h>
31 
32 # include <cmath>
33 
34 DEAL_II_NAMESPACE_OPEN
35 
36 namespace PETScWrappers
37 {
39  {
41  }
42 
43 
44 
46  : solver_control(cn)
48  {}
49 
50 
51 
52  void
54  VectorBase & x,
55  const VectorBase & b,
56  const PreconditionerBase &preconditioner)
57  {
58  /*
59  TODO: PETSc duplicates communicators, so this does not work (you put
60  MPI_COMM_SELF in, but get something other out when you ask PETSc for the
61  communicator. This mainly fails due to the MatrixFree classes, that can not
62  ask PETSc for a communicator. //Timo Heister
63  Assert(A.get_mpi_communicator()==mpi_communicator, ExcMessage("PETSc Solver
64  and Matrix need to use the same MPI_Comm."));
65  Assert(x.get_mpi_communicator()==mpi_communicator, ExcMessage("PETSc Solver
66  and Vector need to use the same MPI_Comm."));
67  Assert(b.get_mpi_communicator()==mpi_communicator, ExcMessage("PETSc Solver
68  and Vector need to use the same MPI_Comm."));
69  */
70 
71  // first create a solver object if this
72  // is necessary
73  if (solver_data.get() == nullptr)
74  {
75  solver_data = std_cxx14::make_unique<SolverData>();
76 
77  PetscErrorCode ierr = KSPCreate(mpi_communicator, &solver_data->ksp);
78  AssertThrow(ierr == 0, ExcPETScError(ierr));
79 
80  // let derived classes set the solver
81  // type, and the preconditioning
82  // object set the type of
83  // preconditioner
85 
86  ierr = KSPSetPC(solver_data->ksp, preconditioner.get_pc());
87  AssertThrow(ierr == 0, ExcPETScError(ierr));
88 
89  // make sure the preconditioner has an associated matrix set
90  const Mat B = preconditioner;
91  AssertThrow(B != nullptr,
92  ExcMessage("PETSc preconditioner should have an"
93  "associated matrix set to be used in solver."));
94 
95  // setting the preconditioner overwrites the used matrices.
96  // hence, we need to set the matrices after the preconditioner.
97 # if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
98  // the last argument is irrelevant here,
99  // since we use the solver only once anyway
100  ierr = KSPSetOperators(solver_data->ksp,
101  A,
102  preconditioner,
103  SAME_PRECONDITIONER);
104 # else
105  ierr = KSPSetOperators(solver_data->ksp, A, preconditioner);
106 # endif
107  AssertThrow(ierr == 0, ExcPETScError(ierr));
108 
109  // then a convergence monitor
110  // function. that function simply
111  // checks with the solver_control
112  // object we have in this object for
113  // convergence
114  ierr = KSPSetConvergenceTest(solver_data->ksp,
116  reinterpret_cast<void *>(&solver_control),
117  nullptr);
118  AssertThrow(ierr == 0, ExcPETScError(ierr));
119  }
120 
121  // set the command line option prefix name
122  PetscErrorCode ierr =
123  KSPSetOptionsPrefix(solver_data->ksp, prefix_name.c_str());
124  AssertThrow(ierr == 0, ExcPETScError(ierr));
125 
126  // set the command line options provided
127  // by the user to override the defaults
128  ierr = KSPSetFromOptions(solver_data->ksp);
129  AssertThrow(ierr == 0, ExcPETScError(ierr));
130 
131  // then do the real work: set up solver
132  // internal data and solve the
133  // system.
134  ierr = KSPSetUp(solver_data->ksp);
135  AssertThrow(ierr == 0, ExcPETScError(ierr));
136 
137  ierr = KSPSolve(solver_data->ksp, b, x);
138  AssertThrow(ierr == 0, ExcPETScError(ierr));
139 
140  // do not destroy solver object
141  // solver_data.reset ();
142 
143  // in case of failure: throw
144  // exception
146  AssertThrow(false,
149  // otherwise exit as normal
150  }
151 
152 
153  void
154  SolverBase::set_prefix(const std::string &prefix)
155  {
156  prefix_name = prefix;
157  }
158 
159 
160  void
162  {
163  solver_data.reset();
164  }
165 
166 
167  SolverControl &
169  {
170  return solver_control;
171  }
172 
173 
174  int
176  const PetscInt iteration,
177  const PetscReal residual_norm,
178  KSPConvergedReason *reason,
179  void * solver_control_x)
180  {
182  *reinterpret_cast<SolverControl *>(solver_control_x);
183 
184  const SolverControl::State state =
185  solver_control.check(iteration, residual_norm);
186 
187  switch (state)
188  {
189  case ::SolverControl::iterate:
190  *reason = KSP_CONVERGED_ITERATING;
191  break;
192 
193  case ::SolverControl::success:
194  *reason = static_cast<KSPConvergedReason>(1);
195  break;
196 
197  case ::SolverControl::failure:
199  *reason = KSP_DIVERGED_ITS;
200  else
201  *reason = KSP_DIVERGED_DTOL;
202  break;
203 
204  default:
205  Assert(false, ExcNotImplemented());
206  }
207 
208  // return without failure
209  return 0;
210  }
211 
212  void
214  {
215  PetscErrorCode ierr;
216 
217  solver_data = std_cxx14::make_unique<SolverData>();
218 
219  ierr = KSPCreate(mpi_communicator, &solver_data->ksp);
220  AssertThrow(ierr == 0, ExcPETScError(ierr));
221 
222  // let derived classes set the solver
223  // type, and the preconditioning
224  // object set the type of
225  // preconditioner
227 
228  ierr = KSPSetPC(solver_data->ksp, preconditioner.get_pc());
229  AssertThrow(ierr == 0, ExcPETScError(ierr));
230 
231  // then a convergence monitor
232  // function. that function simply
233  // checks with the solver_control
234  // object we have in this object for
235  // convergence
236  ierr = KSPSetConvergenceTest(solver_data->ksp,
238  reinterpret_cast<void *>(&solver_control),
239  nullptr);
240  AssertThrow(ierr == 0, ExcPETScError(ierr));
241 
242  // set the command line options provided
243  // by the user to override the defaults
244  ierr = KSPSetFromOptions(solver_data->ksp);
245  AssertThrow(ierr == 0, ExcPETScError(ierr));
246  }
247 
248 
249 
250  /* ---------------------- SolverRichardson ------------------------ */
251 
253  : omega(omega)
254  {}
255 
256 
257 
259  const MPI_Comm & mpi_communicator,
260  const AdditionalData &data)
262  , additional_data(data)
263  {}
264 
265 
266  void
268  {
269  PetscErrorCode ierr = KSPSetType(ksp, KSPRICHARDSON);
270  AssertThrow(ierr == 0, ExcPETScError(ierr));
271 
272  // set the damping factor from the data
273  ierr = KSPRichardsonSetScale(ksp, additional_data.omega);
274  AssertThrow(ierr == 0, ExcPETScError(ierr));
275 
276  // in the deal.II solvers, we always
277  // honor the initial guess in the
278  // solution vector. do so here as well:
279  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
280  AssertThrow(ierr == 0, ExcPETScError(ierr));
281 
282  // Hand over the absolute
283  // tolerance and the maximum
284  // iteration number to the PETSc
285  // convergence criterion. The
286  // custom deal.II SolverControl
287  // object is ignored by the PETSc
288  // Richardson method (when no
289  // PETSc monitoring is present),
290  // since in this case PETSc
291  // uses a faster version of
292  // the Richardson iteration,
293  // where no residual is
294  // available.
295  ierr = KSPSetTolerances(ksp,
296  PETSC_DEFAULT,
297  this->solver_control.tolerance(),
298  PETSC_DEFAULT,
299  this->solver_control.max_steps() + 1);
300  AssertThrow(ierr == 0, ExcPETScError(ierr));
301  }
302 
303 
304  /* ---------------------- SolverChebychev ------------------------ */
305 
307  const MPI_Comm & mpi_communicator,
308  const AdditionalData &data)
309  : SolverBase(cn, mpi_communicator)
310  , additional_data(data)
311  {}
312 
313 
314  void
316  {
317  PetscErrorCode ierr = KSPSetType(ksp, KSPCHEBYSHEV);
318  AssertThrow(ierr == 0, ExcPETScError(ierr));
319 
320  // in the deal.II solvers, we always
321  // honor the initial guess in the
322  // solution vector. do so here as well:
323  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
324  AssertThrow(ierr == 0, ExcPETScError(ierr));
325  }
326 
327 
328  /* ---------------------- SolverCG ------------------------ */
329 
331  const MPI_Comm & mpi_communicator,
332  const AdditionalData &data)
333  : SolverBase(cn, mpi_communicator)
334  , additional_data(data)
335  {}
336 
337 
338  void
340  {
341  PetscErrorCode ierr = KSPSetType(ksp, KSPCG);
342  AssertThrow(ierr == 0, ExcPETScError(ierr));
343 
344  // in the deal.II solvers, we always
345  // honor the initial guess in the
346  // solution vector. do so here as well:
347  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
348  AssertThrow(ierr == 0, ExcPETScError(ierr));
349  }
350 
351 
352  /* ---------------------- SolverBiCG ------------------------ */
353 
355  const MPI_Comm & mpi_communicator,
356  const AdditionalData &data)
357  : SolverBase(cn, mpi_communicator)
358  , additional_data(data)
359  {}
360 
361 
362  void
364  {
365  PetscErrorCode ierr = KSPSetType(ksp, KSPBICG);
366  AssertThrow(ierr == 0, ExcPETScError(ierr));
367 
368  // in the deal.II solvers, we always
369  // honor the initial guess in the
370  // solution vector. do so here as well:
371  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
372  AssertThrow(ierr == 0, ExcPETScError(ierr));
373  }
374 
375 
376  /* ---------------------- SolverGMRES ------------------------ */
377 
379  const unsigned int restart_parameter,
380  const bool right_preconditioning)
381  : restart_parameter(restart_parameter)
382  , right_preconditioning(right_preconditioning)
383  {}
384 
385 
386 
388  const MPI_Comm & mpi_communicator,
389  const AdditionalData &data)
391  , additional_data(data)
392  {}
393 
394 
395  void
397  {
398  PetscErrorCode ierr = KSPSetType(ksp, KSPGMRES);
399  AssertThrow(ierr == 0, ExcPETScError(ierr));
400 
401  // set the restart parameter from the
402  // data. we would like to use the simple
403  // code that is commented out, but this
404  // leads to nasty warning and error
405  // messages due to some stupidity on
406  // PETSc's side: KSPGMRESSetRestart is
407  // implemented as a macro in which return
408  // statements are hidden. This may work
409  // if people strictly follow the PETSc
410  // coding style of always having
411  // functions return an integer error
412  // code, but the present function isn't
413  // like this.
414  /*
415  ierr = KSPGMRESSetRestart (ksp, additional_data.restart_parameter);
416  AssertThrow (ierr == 0, ExcPETScError(ierr));
417  */
418  // so rather expand their macros by hand,
419  // and do some equally nasty stuff that at
420  // least doesn't yield warnings...
421  int (*fun_ptr)(KSP, int);
422  ierr = PetscObjectQueryFunction(reinterpret_cast<PetscObject>(ksp),
423  "KSPGMRESSetRestart_C",
424  reinterpret_cast<void (**)()>(&fun_ptr));
425  AssertThrow(ierr == 0, ExcPETScError(ierr));
426 
427  ierr = (*fun_ptr)(ksp, additional_data.restart_parameter);
428  AssertThrow(ierr == 0, ExcPETScError(ierr));
429 
430  // Set preconditioning side to
431  // right
433  {
434  ierr = KSPSetPCSide(ksp, PC_RIGHT);
435  AssertThrow(ierr == 0, ExcPETScError(ierr));
436  }
437 
438  // in the deal.II solvers, we always
439  // honor the initial guess in the
440  // solution vector. do so here as well:
441  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
442  AssertThrow(ierr == 0, ExcPETScError(ierr));
443  }
444 
445 
446  /* ---------------------- SolverBicgstab ------------------------ */
447 
449  const MPI_Comm & mpi_communicator,
450  const AdditionalData &data)
451  : SolverBase(cn, mpi_communicator)
452  , additional_data(data)
453  {}
454 
455 
456  void
458  {
459  PetscErrorCode ierr = KSPSetType(ksp, KSPBCGS);
460  AssertThrow(ierr == 0, ExcPETScError(ierr));
461 
462  // in the deal.II solvers, we always
463  // honor the initial guess in the
464  // solution vector. do so here as well:
465  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
466  AssertThrow(ierr == 0, ExcPETScError(ierr));
467  }
468 
469 
470  /* ---------------------- SolverCGS ------------------------ */
471 
473  const MPI_Comm & mpi_communicator,
474  const AdditionalData &data)
475  : SolverBase(cn, mpi_communicator)
476  , additional_data(data)
477  {}
478 
479 
480  void
482  {
483  PetscErrorCode ierr = KSPSetType(ksp, KSPCGS);
484  AssertThrow(ierr == 0, ExcPETScError(ierr));
485 
486  // in the deal.II solvers, we always
487  // honor the initial guess in the
488  // solution vector. do so here as well:
489  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
490  AssertThrow(ierr == 0, ExcPETScError(ierr));
491  }
492 
493 
494  /* ---------------------- SolverTFQMR ------------------------ */
495 
497  const MPI_Comm & mpi_communicator,
498  const AdditionalData &data)
499  : SolverBase(cn, mpi_communicator)
500  , additional_data(data)
501  {}
502 
503 
504  void
506  {
507  PetscErrorCode ierr = KSPSetType(ksp, KSPTFQMR);
508  AssertThrow(ierr == 0, ExcPETScError(ierr));
509 
510  // in the deal.II solvers, we always
511  // honor the initial guess in the
512  // solution vector. do so here as well:
513  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
514  AssertThrow(ierr == 0, ExcPETScError(ierr));
515  }
516 
517 
518  /* ---------------------- SolverTCQMR ------------------------ */
519 
521  const MPI_Comm & mpi_communicator,
522  const AdditionalData &data)
523  : SolverBase(cn, mpi_communicator)
524  , additional_data(data)
525  {}
526 
527 
528  void
530  {
531  PetscErrorCode ierr = KSPSetType(ksp, KSPTCQMR);
532  AssertThrow(ierr == 0, ExcPETScError(ierr));
533 
534  // in the deal.II solvers, we always
535  // honor the initial guess in the
536  // solution vector. do so here as well:
537  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
538  AssertThrow(ierr == 0, ExcPETScError(ierr));
539  }
540 
541 
542  /* ---------------------- SolverCR ------------------------ */
543 
545  const MPI_Comm & mpi_communicator,
546  const AdditionalData &data)
547  : SolverBase(cn, mpi_communicator)
548  , additional_data(data)
549  {}
550 
551 
552  void
554  {
555  PetscErrorCode ierr = KSPSetType(ksp, KSPCR);
556  AssertThrow(ierr == 0, ExcPETScError(ierr));
557 
558  // in the deal.II solvers, we always
559  // honor the initial guess in the
560  // solution vector. do so here as well:
561  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
562  AssertThrow(ierr == 0, ExcPETScError(ierr));
563  }
564 
565 
566  /* ---------------------- SolverLSQR ------------------------ */
567 
569  const MPI_Comm & mpi_communicator,
570  const AdditionalData &data)
571  : SolverBase(cn, mpi_communicator)
572  , additional_data(data)
573  {}
574 
575 
576  void
578  {
579  PetscErrorCode ierr = KSPSetType(ksp, KSPLSQR);
580  AssertThrow(ierr == 0, ExcPETScError(ierr));
581 
582  // in the deal.II solvers, we always
583  // honor the initial guess in the
584  // solution vector. do so here as well:
585  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_TRUE);
586  AssertThrow(ierr == 0, ExcPETScError(ierr));
587  }
588 
589 
590  /* ---------------------- SolverPreOnly ------------------------ */
591 
593  const MPI_Comm & mpi_communicator,
594  const AdditionalData &data)
595  : SolverBase(cn, mpi_communicator)
596  , additional_data(data)
597  {}
598 
599 
600  void
602  {
603  PetscErrorCode ierr = KSPSetType(ksp, KSPPREONLY);
604  AssertThrow(ierr == 0, ExcPETScError(ierr));
605 
606  // The KSPPREONLY solver of
607  // PETSc never calls the convergence
608  // monitor, which leads to failure
609  // even when everything was ok.
610  // Therefore the SolverControl status
611  // is set to some nice values, which
612  // guarantee a nice result at the end
613  // of the solution process.
614  solver_control.check(1, 0.0);
615 
616  // Using the PREONLY solver with
617  // a nonzero initial guess leads
618  // PETSc to produce some error messages.
619  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_FALSE);
620  AssertThrow(ierr == 0, ExcPETScError(ierr));
621  }
622 
623 
624  /* ---------------------- SparseDirectMUMPS------------------------ */
625 
627  {
629  // the 'pc' object is owned by the 'ksp' object, and consequently
630  // does not have to be destroyed explicitly here
631  }
632 
633 
635  const MPI_Comm & mpi_communicator,
636  const AdditionalData &data)
638  , additional_data(data)
639  , symmetric_mode(false)
640  {}
641 
642 
643  void
645  {
651  PetscErrorCode ierr = KSPSetType(ksp, KSPPREONLY);
652  AssertThrow(ierr == 0, ExcPETScError(ierr));
653 
660  solver_control.check(1, 0.0);
661 
666  ierr = KSPSetInitialGuessNonzero(ksp, PETSC_FALSE);
667  AssertThrow(ierr == 0, ExcPETScError(ierr));
668  }
669 
670  void
672  VectorBase & x,
673  const VectorBase &b)
674  {
675 # ifdef DEAL_II_PETSC_WITH_MUMPS
676 
679  Mat F;
680 
686  PetscInt ival = 2, icntl = 7;
690  PetscInt its;
694  PetscReal rnorm;
695 
699  if (solver_data == nullptr)
700  {
701  solver_data = std_cxx14::make_unique<SolverDataMUMPS>();
702 
707  PetscErrorCode ierr = KSPCreate(mpi_communicator, &solver_data->ksp);
708  AssertThrow(ierr == 0, ExcPETScError(ierr));
709 
714 # if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
715  ierr =
716  KSPSetOperators(solver_data->ksp, A, A, DIFFERENT_NONZERO_PATTERN);
717 # else
718  ierr = KSPSetOperators(solver_data->ksp, A, A);
719 # endif
720  AssertThrow(ierr == 0, ExcPETScError(ierr));
721 
725  set_solver_type(solver_data->ksp);
726 
730  ierr = KSPGetPC(solver_data->ksp, &solver_data->pc);
731  AssertThrow(ierr == 0, ExcPETScError(ierr));
732 
737  if (symmetric_mode)
738  ierr = PCSetType(solver_data->pc, PCCHOLESKY);
739  else
740  ierr = PCSetType(solver_data->pc, PCLU);
741  AssertThrow(ierr == 0, ExcPETScError(ierr));
742 
747  ierr = KSPSetConvergenceTest(solver_data->ksp,
749  reinterpret_cast<void *>(&solver_control),
750  PETSC_NULL);
751  AssertThrow(ierr == 0, ExcPETScError(ierr));
752 
758 # if DEAL_II_PETSC_VERSION_LT(3, 9, 0)
759  ierr = PCFactorSetMatSolverPackage(solver_data->pc, MATSOLVERMUMPS);
760 # else
761  ierr = PCFactorSetMatSolverType(solver_data->pc, MATSOLVERMUMPS);
762 # endif
763  AssertThrow(ierr == 0, ExcPETScError(ierr));
764 
768 # if DEAL_II_PETSC_VERSION_LT(3, 9, 0)
769  ierr = PCFactorSetUpMatSolverPackage(solver_data->pc);
770 # else
771  ierr = PCFactorSetUpMatSolverType(solver_data->pc);
772 # endif
773  AssertThrow(ierr == 0, ExcPETScError(ierr));
774 
780  ierr = PCFactorGetMatrix(solver_data->pc, &F);
781  AssertThrow(ierr == 0, ExcPETScError(ierr));
782 
786  ierr = MatMumpsSetIcntl(F, icntl, ival);
787  AssertThrow(ierr == 0, ExcPETScError(ierr));
788 
792  ierr = KSPSetOptionsPrefix(solver_data->ksp, prefix_name.c_str());
793  AssertThrow(ierr == 0, ExcPETScError(ierr));
794 
799  ierr = KSPSetFromOptions(solver_data->ksp);
800  AssertThrow(ierr == 0, ExcPETScError(ierr));
801  }
802 
806  PetscErrorCode ierr = KSPSolve(solver_data->ksp, b, x);
807  AssertThrow(ierr == 0, ExcPETScError(ierr));
808 
813  {
814  AssertThrow(false,
817  }
818  else
819  {
824  ierr = KSPGetIterationNumber(solver_data->ksp, &its);
825  AssertThrow(ierr == 0, ExcPETScError(ierr));
826  ierr = KSPGetResidualNorm(solver_data->ksp, &rnorm);
827  AssertThrow(ierr == 0, ExcPETScError(ierr));
828  }
829 
830 # else // DEAL_II_PETSC_WITH_MUMPS
831  Assert(
832  false,
833  ExcMessage(
834  "Your PETSc installation does not include a copy of "
835  "the MUMPS package necessary for this solver. You will need to configure "
836  "PETSc so that it includes MUMPS, recompile it, and then re-configure "
837  "and recompile deal.II as well."));
838 
839  // Cast to void to silence compiler warnings
840  (void)A;
841  (void)x;
842  (void)b;
843 # endif
844  }
845 
846  PetscErrorCode
848  const PetscInt iteration,
849  const PetscReal residual_norm,
850  KSPConvergedReason *reason,
851  void * solver_control_x)
852  {
854  *reinterpret_cast<SolverControl *>(solver_control_x);
855 
856  const SolverControl::State state =
857  solver_control.check(iteration, residual_norm);
858 
859  switch (state)
860  {
861  case ::SolverControl::iterate:
862  *reason = KSP_CONVERGED_ITERATING;
863  break;
864 
865  case ::SolverControl::success:
866  *reason = static_cast<KSPConvergedReason>(1);
867  break;
868 
869  case ::SolverControl::failure:
871  *reason = KSP_DIVERGED_ITS;
872  else
873  *reason = KSP_DIVERGED_DTOL;
874  break;
875 
876  default:
877  Assert(false, ExcNotImplemented());
878  }
879 
880  return 0;
881  }
882 
883  void
885  {
886  symmetric_mode = flag;
887  }
888 
889 } // namespace PETScWrappers
890 
891 DEAL_II_NAMESPACE_CLOSE
892 
893 #endif // DEAL_II_WITH_PETSC
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const =0
virtual void set_solver_type(KSP &ksp) const override
virtual State check(const unsigned int step, const double check_value)
virtual void set_solver_type(KSP &ksp) const override
const AdditionalData additional_data
Definition: petsc_solver.h:310
void set_prefix(const std::string &prefix)
AdditionalData(const unsigned int restart_parameter=30, const bool right_preconditioning=false)
SolverRichardson(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:532
void initialize(const PreconditionerBase &preconditioner)
#define AssertThrow(cond, exc)
Definition: exceptions.h:1519
virtual void set_solver_type(KSP &ksp) const override
SparseDirectMUMPS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static ::ExceptionBase & ExcMessage(std::string arg1)
const MPI_Comm mpi_communicator
Definition: petsc_solver.h:182
Stop iteration, goal reached.
std::unique_ptr< SolverData > solver_data
Definition: petsc_solver.h:247
virtual void set_solver_type(KSP &ksp) const override
#define Assert(cond, exc)
Definition: exceptions.h:1407
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
Definition: petsc_solver.cc:53
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
Definition: petsc_solver.cc:45
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override
SolverCR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverControl & solver_control
Definition: petsc_solver.h:177
double last_value() const
SolverBiCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverLSQR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverChebychev(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
double tolerance() const
SolverCGS(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
unsigned int last_step() const
PetscErrorCode destroy_krylov_solver(KSP &krylov_solver)
SolverBicgstab(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
const AdditionalData additional_data
Definition: petsc_solver.h:967
SolverControl & control() const
virtual void set_solver_type(KSP &ksp) const override
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b)
SolverTCQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
State last_check() const
unsigned int max_steps() const
SolverTFQMR(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
static ::ExceptionBase & ExcNotImplemented()
void set_symmetric_mode(const bool flag)
static PetscErrorCode convergence_test(KSP ksp, const PetscInt iteration, const PetscReal residual_norm, KSPConvergedReason *reason, void *solver_control)
SolverCG(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverPreOnly(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
SolverGMRES(SolverControl &cn, const MPI_Comm &mpi_communicator=PETSC_COMM_SELF, const AdditionalData &data=AdditionalData())
virtual void set_solver_type(KSP &ksp) const override
virtual void set_solver_type(KSP &ksp) const override