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