Reference documentation for deal.II version 9.0.0
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 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/base/std_cxx14/memory.h>
19 #include <deal.II/lac/petsc_solver.h>
20 
21 #ifdef DEAL_II_WITH_PETSC
22 
23 # include <deal.II/lac/exceptions.h>
24 # include <deal.II/lac/petsc_compatibility.h>
25 # include <deal.II/lac/petsc_matrix_base.h>
26 # include <deal.II/lac/petsc_vector_base.h>
27 # include <deal.II/lac/petsc_precondition.h>
28 # include <cmath>
29 
30 #include <petscversion.h>
31 
32 DEAL_II_NAMESPACE_OPEN
33 
34 namespace PETScWrappers
35 {
36 
38  {
40  }
41 
42 
43 
45  const MPI_Comm &mpi_communicator)
46  :
47  solver_control (cn),
49  {}
50 
51 
52 
53  void
55  VectorBase &x,
56  const VectorBase &b,
57  const PreconditionerBase &preconditioner)
58  {
59  /*
60  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
61  Assert(A.get_mpi_communicator()==mpi_communicator, ExcMessage("PETSc Solver and Matrix need to use the same MPI_Comm."));
62  Assert(x.get_mpi_communicator()==mpi_communicator, ExcMessage("PETSc Solver and Vector need to use the same MPI_Comm."));
63  Assert(b.get_mpi_communicator()==mpi_communicator, ExcMessage("PETSc Solver and Vector need to use the same MPI_Comm."));
64  */
65 
66  // first create a solver object if this
67  // is necessary
68  if (solver_data.get() == nullptr)
69  {
70  solver_data = std_cxx14::make_unique<SolverData>();
71 
72  PetscErrorCode ierr = KSPCreate (mpi_communicator, &solver_data->ksp);
73  AssertThrow (ierr == 0, ExcPETScError(ierr));
74 
75  // let derived classes set the solver
76  // type, and the preconditioning
77  // object set the type of
78  // preconditioner
80 
81  ierr = KSPSetPC (solver_data->ksp, preconditioner.get_pc());
82  AssertThrow (ierr == 0, ExcPETScError(ierr));
83 
84  // make sure the preconditioner has an associated matrix set
85  const Mat B = preconditioner;
86  AssertThrow (B != nullptr,
87  ExcMessage("PETSc preconditioner should have an"
88  "associated matrix set to be used in solver."));
89 
90  // setting the preconditioner overwrites the used matrices.
91  // hence, we need to set the matrices after the preconditioner.
92 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
93  // the last argument is irrelevant here,
94  // since we use the solver only once anyway
95  ierr = KSPSetOperators (solver_data->ksp, A, preconditioner,
96  SAME_PRECONDITIONER);
97 #else
98  ierr = KSPSetOperators (solver_data->ksp, A, preconditioner);
99 #endif
100  AssertThrow (ierr == 0, ExcPETScError(ierr));
101 
102  // then a convergence monitor
103  // function. that function simply
104  // checks with the solver_control
105  // object we have in this object for
106  // convergence
107  ierr = KSPSetConvergenceTest (solver_data->ksp, &convergence_test,
108  reinterpret_cast<void *>(&solver_control),
109  nullptr);
110  AssertThrow (ierr == 0, ExcPETScError(ierr));
111  }
112 
113  // set the command line option prefix name
114  PetscErrorCode ierr = KSPSetOptionsPrefix(solver_data->ksp, prefix_name.c_str());
115  AssertThrow (ierr == 0, ExcPETScError(ierr));
116 
117  // set the command line options provided
118  // by the user to override the defaults
119  ierr = KSPSetFromOptions (solver_data->ksp);
120  AssertThrow (ierr == 0, ExcPETScError(ierr));
121 
122  // then do the real work: set up solver
123  // internal data and solve the
124  // system.
125  ierr = KSPSetUp (solver_data->ksp);
126  AssertThrow (ierr == 0, ExcPETScError(ierr));
127 
128  ierr = KSPSolve (solver_data->ksp, b, x);
129  AssertThrow (ierr == 0, ExcPETScError(ierr));
130 
131  // do not destroy solver object
132 // solver_data.reset ();
133 
134  // in case of failure: throw
135  // exception
139  // otherwise exit as normal
140  }
141 
142 
143  void
144  SolverBase::set_prefix(const std::string &prefix)
145  {
146  prefix_name = prefix ;
147  }
148 
149 
150  void
152  {
153  solver_data.reset ();
154  }
155 
156 
157  SolverControl &
159  {
160  return solver_control;
161  }
162 
163 
164  int
166  const PetscInt iteration,
167  const PetscReal residual_norm,
168  KSPConvergedReason *reason,
169  void *solver_control_x)
170  {
171  SolverControl &solver_control = *reinterpret_cast<SolverControl *>(solver_control_x);
172 
173  const SolverControl::State state
174  = solver_control.check (iteration, residual_norm);
175 
176  switch (state)
177  {
178  case ::SolverControl::iterate:
179  *reason = KSP_CONVERGED_ITERATING;
180  break;
181 
182  case ::SolverControl::success:
183  *reason = static_cast<KSPConvergedReason>(1);
184  break;
185 
186  case ::SolverControl::failure:
188  *reason = KSP_DIVERGED_ITS;
189  else
190  *reason = KSP_DIVERGED_DTOL;
191  break;
192 
193  default:
194  Assert (false, ExcNotImplemented());
195  }
196 
197  // return without failure
198  return 0;
199  }
200 
201  void
203  {
204  PetscErrorCode ierr;
205 
206  solver_data = std_cxx14::make_unique<SolverData>();
207 
208  ierr = KSPCreate (mpi_communicator, &solver_data->ksp);
209  AssertThrow (ierr == 0, ExcPETScError(ierr));
210 
211  // let derived classes set the solver
212  // type, and the preconditioning
213  // object set the type of
214  // preconditioner
216 
217  ierr = KSPSetPC (solver_data->ksp, preconditioner.get_pc());
218  AssertThrow (ierr == 0, ExcPETScError(ierr));
219 
220  // then a convergence monitor
221  // function. that function simply
222  // checks with the solver_control
223  // object we have in this object for
224  // convergence
225  ierr = KSPSetConvergenceTest (solver_data->ksp, &convergence_test,
226  reinterpret_cast<void *>(&solver_control),
227  nullptr);
228  AssertThrow (ierr == 0, ExcPETScError(ierr));
229 
230  // set the command line options provided
231  // by the user to override the defaults
232  ierr = KSPSetFromOptions (solver_data->ksp);
233  AssertThrow (ierr == 0, ExcPETScError(ierr));
234  }
235 
236 
237 
238  /* ---------------------- SolverRichardson ------------------------ */
239 
241  AdditionalData (const double omega)
242  :
243  omega (omega)
244  {}
245 
246 
247 
249  const MPI_Comm &mpi_communicator,
250  const AdditionalData &data)
251  :
253  additional_data (data)
254  {}
255 
256 
257  void
259  {
260  PetscErrorCode ierr = KSPSetType (ksp, KSPRICHARDSON);
261  AssertThrow (ierr == 0, ExcPETScError(ierr));
262 
263  // set the damping factor from the data
264  ierr = KSPRichardsonSetScale (ksp, additional_data.omega);
265  AssertThrow (ierr == 0, ExcPETScError(ierr));
266 
267  // in the deal.II solvers, we always
268  // honor the initial guess in the
269  // solution vector. do so here as well:
270  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
271  AssertThrow (ierr == 0, ExcPETScError(ierr));
272 
273  // Hand over the absolute
274  // tolerance and the maximum
275  // iteration number to the PETSc
276  // convergence criterion. The
277  // custom deal.II SolverControl
278  // object is ignored by the PETSc
279  // Richardson method (when no
280  // PETSc monitoring is present),
281  // since in this case PETSc
282  // uses a faster version of
283  // the Richardson iteration,
284  // where no residual is
285  // available.
286  ierr = KSPSetTolerances(ksp, PETSC_DEFAULT, this->solver_control.tolerance(),
287  PETSC_DEFAULT, this->solver_control.max_steps()+1);
288  AssertThrow (ierr == 0, ExcPETScError(ierr));
289  }
290 
291 
292  /* ---------------------- SolverChebychev ------------------------ */
293 
295  const MPI_Comm &mpi_communicator,
296  const AdditionalData &data)
297  :
298  SolverBase (cn, mpi_communicator),
299  additional_data (data)
300  {}
301 
302 
303  void
305  {
306  PetscErrorCode ierr = KSPSetType (ksp, KSPCHEBYSHEV);
307  AssertThrow (ierr == 0, ExcPETScError(ierr));
308 
309  // in the deal.II solvers, we always
310  // honor the initial guess in the
311  // solution vector. do so here as well:
312  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
313  AssertThrow (ierr == 0, ExcPETScError(ierr));
314  }
315 
316 
317  /* ---------------------- SolverCG ------------------------ */
318 
320  const MPI_Comm &mpi_communicator,
321  const AdditionalData &data)
322  :
323  SolverBase (cn, mpi_communicator),
324  additional_data (data)
325  {}
326 
327 
328  void
329  SolverCG::set_solver_type (KSP &ksp) const
330  {
331  PetscErrorCode ierr = KSPSetType (ksp, KSPCG);
332  AssertThrow (ierr == 0, ExcPETScError(ierr));
333 
334  // in the deal.II solvers, we always
335  // honor the initial guess in the
336  // solution vector. do so here as well:
337  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
338  AssertThrow (ierr == 0, ExcPETScError(ierr));
339  }
340 
341 
342  /* ---------------------- SolverBiCG ------------------------ */
343 
345  const MPI_Comm &mpi_communicator,
346  const AdditionalData &data)
347  :
348  SolverBase (cn, mpi_communicator),
349  additional_data (data)
350  {}
351 
352 
353  void
355  {
356  PetscErrorCode ierr = KSPSetType (ksp, KSPBICG);
357  AssertThrow (ierr == 0, ExcPETScError(ierr));
358 
359  // in the deal.II solvers, we always
360  // honor the initial guess in the
361  // solution vector. do so here as well:
362  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
363  AssertThrow (ierr == 0, ExcPETScError(ierr));
364  }
365 
366 
367  /* ---------------------- SolverGMRES ------------------------ */
368 
370  AdditionalData (const unsigned int restart_parameter,
371  const bool right_preconditioning)
372  :
373  restart_parameter (restart_parameter),
374  right_preconditioning (right_preconditioning)
375  {}
376 
377 
378 
380  const MPI_Comm &mpi_communicator,
381  const AdditionalData &data)
382  :
384  additional_data (data)
385  {}
386 
387 
388  void
390  {
391  PetscErrorCode ierr = KSPSetType (ksp, KSPGMRES);
392  AssertThrow (ierr == 0, ExcPETScError(ierr));
393 
394  // set the restart parameter from the
395  // data. we would like to use the simple
396  // code that is commented out, but this
397  // leads to nasty warning and error
398  // messages due to some stupidity on
399  // PETSc's side: KSPGMRESSetRestart is
400  // implemented as a macro in which return
401  // statements are hidden. This may work
402  // if people strictly follow the PETSc
403  // coding style of always having
404  // functions return an integer error
405  // code, but the present function isn't
406  // like this.
407  /*
408  ierr = KSPGMRESSetRestart (ksp, additional_data.restart_parameter);
409  AssertThrow (ierr == 0, ExcPETScError(ierr));
410  */
411  // so rather expand their macros by hand,
412  // and do some equally nasty stuff that at
413  // least doesn't yield warnings...
414  int (*fun_ptr)(KSP,int);
415  ierr = PetscObjectQueryFunction((PetscObject)(ksp),
416  "KSPGMRESSetRestart_C",
417  (void (* *)())&fun_ptr);
418  AssertThrow (ierr == 0, ExcPETScError(ierr));
419 
420  ierr = (*fun_ptr)(ksp,additional_data.restart_parameter);
421  AssertThrow (ierr == 0, ExcPETScError(ierr));
422 
423  // Set preconditioning side to
424  // right
426  {
427  ierr = KSPSetPCSide(ksp, PC_RIGHT);
428  AssertThrow (ierr == 0, ExcPETScError(ierr));
429  }
430 
431  // in the deal.II solvers, we always
432  // honor the initial guess in the
433  // solution vector. do so here as well:
434  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
435  AssertThrow (ierr == 0, ExcPETScError(ierr));
436  }
437 
438 
439  /* ---------------------- SolverBicgstab ------------------------ */
440 
442  const MPI_Comm &mpi_communicator,
443  const AdditionalData &data)
444  :
445  SolverBase (cn, mpi_communicator),
446  additional_data (data)
447  {}
448 
449 
450  void
452  {
453  PetscErrorCode ierr = KSPSetType (ksp, KSPBCGS);
454  AssertThrow (ierr == 0, ExcPETScError(ierr));
455 
456  // in the deal.II solvers, we always
457  // honor the initial guess in the
458  // solution vector. do so here as well:
459  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
460  AssertThrow (ierr == 0, ExcPETScError(ierr));
461  }
462 
463 
464  /* ---------------------- SolverCGS ------------------------ */
465 
467  const MPI_Comm &mpi_communicator,
468  const AdditionalData &data)
469  :
470  SolverBase (cn, mpi_communicator),
471  additional_data (data)
472  {}
473 
474 
475  void
476  SolverCGS::set_solver_type (KSP &ksp) const
477  {
478  PetscErrorCode ierr = KSPSetType (ksp, KSPCGS);
479  AssertThrow (ierr == 0, ExcPETScError(ierr));
480 
481  // in the deal.II solvers, we always
482  // honor the initial guess in the
483  // solution vector. do so here as well:
484  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
485  AssertThrow (ierr == 0, ExcPETScError(ierr));
486  }
487 
488 
489  /* ---------------------- SolverTFQMR ------------------------ */
490 
492  const MPI_Comm &mpi_communicator,
493  const AdditionalData &data)
494  :
495  SolverBase (cn, mpi_communicator),
496  additional_data (data)
497  {}
498 
499 
500  void
502  {
503  PetscErrorCode ierr = KSPSetType (ksp, KSPTFQMR);
504  AssertThrow (ierr == 0, ExcPETScError(ierr));
505 
506  // in the deal.II solvers, we always
507  // honor the initial guess in the
508  // solution vector. do so here as well:
509  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
510  AssertThrow (ierr == 0, ExcPETScError(ierr));
511  }
512 
513 
514  /* ---------------------- SolverTCQMR ------------------------ */
515 
517  const MPI_Comm &mpi_communicator,
518  const AdditionalData &data)
519  :
520  SolverBase (cn, mpi_communicator),
521  additional_data (data)
522  {}
523 
524 
525  void
527  {
528  PetscErrorCode ierr = KSPSetType (ksp, KSPTCQMR);
529  AssertThrow (ierr == 0, ExcPETScError(ierr));
530 
531  // in the deal.II solvers, we always
532  // honor the initial guess in the
533  // solution vector. do so here as well:
534  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
535  AssertThrow (ierr == 0, ExcPETScError(ierr));
536  }
537 
538 
539  /* ---------------------- SolverCR ------------------------ */
540 
542  const MPI_Comm &mpi_communicator,
543  const AdditionalData &data)
544  :
545  SolverBase (cn, mpi_communicator),
546  additional_data (data)
547  {}
548 
549 
550  void
551  SolverCR::set_solver_type (KSP &ksp) const
552  {
553  PetscErrorCode ierr = KSPSetType (ksp, KSPCR);
554  AssertThrow (ierr == 0, ExcPETScError(ierr));
555 
556  // in the deal.II solvers, we always
557  // honor the initial guess in the
558  // solution vector. do so here as well:
559  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
560  AssertThrow (ierr == 0, ExcPETScError(ierr));
561  }
562 
563 
564  /* ---------------------- SolverLSQR ------------------------ */
565 
567  const MPI_Comm &mpi_communicator,
568  const AdditionalData &data)
569  :
570  SolverBase (cn, mpi_communicator),
571  additional_data (data)
572  {}
573 
574 
575  void
577  {
578  PetscErrorCode ierr = KSPSetType (ksp, KSPLSQR);
579  AssertThrow (ierr == 0, ExcPETScError(ierr));
580 
581  // in the deal.II solvers, we always
582  // honor the initial guess in the
583  // solution vector. do so here as well:
584  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_TRUE);
585  AssertThrow (ierr == 0, ExcPETScError(ierr));
586  }
587 
588 
589  /* ---------------------- SolverPreOnly ------------------------ */
590 
592  const MPI_Comm &mpi_communicator,
593  const AdditionalData &data)
594  :
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  {
628  destroy_krylov_solver (ksp);
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)
637  :
639  additional_data (data),
640  symmetric_mode(false)
641  {}
642 
643 
644  void
646  {
652  PetscErrorCode ierr = KSPSetType (ksp, KSPPREONLY);
653  AssertThrow (ierr == 0, ExcPETScError(ierr));
654 
661  solver_control.check (1, 0.0);
662 
667  ierr = KSPSetInitialGuessNonzero (ksp, PETSC_FALSE);
668  AssertThrow (ierr == 0, ExcPETScError(ierr));
669  }
670 
671  void
673  VectorBase &x,
674  const VectorBase &b)
675  {
676 #ifdef DEAL_II_PETSC_WITH_MUMPS
677 
680  Mat F;
681 
687  PetscInt ival=2, icntl=7;
691  PetscInt its;
695  PetscReal rnorm;
696 
700  if (solver_data == nullptr)
701  {
702  solver_data = std_cxx14::make_unique<SolverDataMUMPS >();
703 
708  PetscErrorCode ierr = KSPCreate (mpi_communicator, &solver_data->ksp);
709  AssertThrow (ierr == 0, ExcPETScError(ierr));
710 
715 #if DEAL_II_PETSC_VERSION_LT(3, 5, 0)
716  ierr = KSPSetOperators (solver_data->ksp, A, A,
717  DIFFERENT_NONZERO_PATTERN);
718 #else
719  ierr = KSPSetOperators (solver_data->ksp, A, A);
720 #endif
721  AssertThrow (ierr == 0, ExcPETScError(ierr));
722 
726  set_solver_type (solver_data->ksp);
727 
731  ierr = KSPGetPC (solver_data->ksp, & solver_data->pc);
732  AssertThrow (ierr == 0, ExcPETScError(ierr));
733 
738  if (symmetric_mode)
739  ierr = PCSetType (solver_data->pc, PCCHOLESKY);
740  else
741  ierr = PCSetType (solver_data->pc, PCLU);
742  AssertThrow (ierr == 0, ExcPETScError(ierr));
743 
748  ierr = KSPSetConvergenceTest (solver_data->ksp, &convergence_test,
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  }
803 
807  PetscErrorCode ierr = KSPSolve (solver_data->ksp, b, x);
808  AssertThrow (ierr == 0, ExcPETScError(ierr));
809 
814  {
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 (false,
832  ExcMessage ("Your PETSc installation does not include a copy of "
833  "the MUMPS package necessary for this solver. You will need to configure "
834  "PETSc so that it includes MUMPS, recompile it, and then re-configure "
835  "and recompile deal.II as well."));
836 
837  // Cast to void to silence compiler warnings
838  (void) A;
839  (void) x;
840  (void) b;
841 #endif
842 
843  }
844 
845  PetscErrorCode SparseDirectMUMPS::convergence_test (KSP /*ksp*/,
846  const PetscInt iteration,
847  const PetscReal residual_norm,
848  KSPConvergedReason *reason,
849  void *solver_control_x)
850  {
851  SolverControl &solver_control = *reinterpret_cast<SolverControl *>(solver_control_x);
852 
853  const SolverControl::State state
854  = solver_control.check (iteration, residual_norm);
855 
856  switch (state)
857  {
858  case ::SolverControl::iterate:
859  *reason = KSP_CONVERGED_ITERATING;
860  break;
861 
862  case ::SolverControl::success:
863  *reason = static_cast<KSPConvergedReason>(1);
864  break;
865 
866  case ::SolverControl::failure:
868  *reason = KSP_DIVERGED_ITS;
869  else
870  *reason = KSP_DIVERGED_DTOL;
871  break;
872 
873  default:
874  Assert (false, ExcNotImplemented());
875  }
876 
877  return 0;
878  }
879 
880  void
882  {
883  symmetric_mode = flag;
884  }
885 
886 }
887 
888 DEAL_II_NAMESPACE_CLOSE
889 
890 #endif // DEAL_II_WITH_PETSC
void solve(const MatrixBase &A, VectorBase &x, const VectorBase &b, const PreconditionerBase &preconditioner)
Definition: petsc_solver.cc:54
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:307
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
const AdditionalData additional_data
Definition: petsc_solver.h:525
virtual void set_solver_type(KSP &ksp) const
#define AssertThrow(cond, exc)
Definition: exceptions.h:1221
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:179
Stop iteration, goal reached.
std::unique_ptr< SolverData > solver_data
Definition: petsc_solver.h:243
#define Assert(cond, exc)
Definition: exceptions.h:1142
SolverBase(SolverControl &cn, const MPI_Comm &mpi_communicator)
Definition: petsc_solver.cc:44
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:174
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:953
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