www.mooseframework.org
Classes | Functions
Moose::PetscSupport Namespace Reference

Classes

class  PetscOptions
 A struct for storing the various types of petsc options and values. More...
 

Functions

void petscSetOptions (const PetscOptions &po, const SolverParams &solver_params)
 A function for setting the PETSc options in PETSc from the options supplied to MOOSE. More...
 
void petscSetKSPDefaults (FEProblemBase &problem, KSP ksp)
 Set the default options for a KSP. More...
 
template<typename T >
void setLinearSolverDefaults (FEProblemBase &problem, LinearSolver< T > &linear_solver)
 Set the defaults for a libMesh LinearSolver. More...
 
void petscSetDefaults (FEProblemBase &problem)
 Sets the default options for PETSc. More...
 
void petscSetupDM (NonlinearSystemBase &nl, const std::string &dm_name)
 Setup the PETSc DM object. More...
 
PetscErrorCode petscSetupOutput (CommandLine *cmd_line)
 
void outputNorm (libMesh::Real old_norm, libMesh::Real norm, bool use_color=false)
 Helper function for outputting the norm values with/without color. More...
 
PetscErrorCode petscLinearMonitor (KSP, PetscInt its, PetscReal rnorm, void *void_ptr)
 Helper function for displaying the linear residual during PETSC solve. More...
 
void storePetscOptions (FEProblemBase &fe_problem, const InputParameters &params)
 Stores the PETSc options supplied from the InputParameters with MOOSE. More...
 
void processPetscFlags (const MultiMooseEnum &petsc_flags, PetscOptions &petsc_options)
 Populate flags in a given PetscOptions object using a vector of input arguments. More...
 
void processPetscPairs (const std::vector< std::pair< MooseEnumItem, std::string >> &petsc_pair_options, const unsigned int mesh_dimension, PetscOptions &petsc_options)
 Populate name and value pairs in a given PetscOptions object using vectors of input arguments. More...
 
std::set< std::string > getPetscValidLineSearches ()
 Returns the valid petsc line search options as a set of strings. More...
 
InputParameters getPetscValidParams ()
 Returns the PETSc options that are common between Executioners and Preconditioners. More...
 
MultiMooseEnum getCommonPetscFlags ()
 A helper function to produce a MultiMooseEnum with commonly used PETSc single options (flags) More...
 
MultiMooseEnum getCommonPetscKeys ()
 A helper function to produce a MultiMooseEnum with commonly used PETSc iname options (keys in key-value pairs) More...
 
bool isSNESVI (FEProblemBase &fe_problem)
 check if SNES type is variational inequalities (VI) solver More...
 
void setSinglePetscOption (const std::string &name, const std::string &value="")
 A wrapper function for dealing with different versions of PetscOptionsSetValue. More...
 
void addPetscOptionsFromCommandline ()
 
void petscSetDefaultPCSide (FEProblemBase &problem, KSP ksp)
 Setup which side we want to apply preconditioner. More...
 
void petscSetDefaultKSPNormType (FEProblemBase &problem, KSP ksp)
 Set norm type. More...
 
void colorAdjacencyMatrix (PetscScalar *adjacency_matrix, unsigned int size, unsigned int colors, std::vector< unsigned int > &vertex_colors, const char *coloring_algorithm)
 This method takes an adjacency matrix, and a desired number of colors and applies a graph coloring algorithm to produce a coloring. More...
 
void disableNonlinearConvergedReason (FEProblemBase &fe_problem)
 disable printing of the nonlinear convergence reason More...
 
void disableLinearConvergedReason (FEProblemBase &fe_problem)
 disable printing of the linear convergence reason More...
 
std::string stringify (const LineSearchType &t)
 
std::string stringify (const MffdType &t)
 
void setSolverOptions (const SolverParams &solver_params)
 
PetscErrorCode petscNonlinearConverged (SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
 
PCSide getPetscPCSide (Moose::PCSideType pcs)
 
KSPNormType getPetscKSPNormType (Moose::MooseKSPNormType kspnorm)
 

Function Documentation

◆ addPetscOptionsFromCommandline()

void Moose::PetscSupport::addPetscOptionsFromCommandline ( )

Definition at line 221 of file PetscSupport.C.

Referenced by petscSetOptions(), and Moose::SlepcSupport::slepcSetOptions().

222 {
223  // commandline options always win
224  // the options from a user commandline will overwrite the existing ones if any conflicts
225  { // Get any options specified on the command-line
226  int argc;
227  char ** args;
228 
229  PetscGetArgs(&argc, &args);
230 #if PETSC_VERSION_LESS_THAN(3, 7, 0)
231  PetscOptionsInsert(&argc, &args, NULL);
232 #else
233  PetscOptionsInsert(LIBMESH_PETSC_NULLPTR, &argc, &args, NULL);
234 #endif
235  }
236 }

◆ colorAdjacencyMatrix()

void Moose::PetscSupport::colorAdjacencyMatrix ( PetscScalar *  adjacency_matrix,
unsigned int  size,
unsigned int  colors,
std::vector< unsigned int > &  vertex_colors,
const char *  coloring_algorithm 
)

This method takes an adjacency matrix, and a desired number of colors and applies a graph coloring algorithm to produce a coloring.

The coloring is returned as a vector of unsigned integers indicating which color or group each vextex in the adjacency matrix belongs to.

Definition at line 953 of file PetscSupport.C.

958 {
959  // Mat A will be a dense matrix from the incoming data structure
960  Mat A;
961  MatCreate(MPI_COMM_SELF, &A);
962  MatSetSizes(A, size, size, size, size);
963  MatSetType(A, MATSEQDENSE);
964  // PETSc requires a non-const data array to populate the matrix
965  MatSeqDenseSetPreallocation(A, adjacency_matrix);
966  MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
967  MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
968 
969  // Convert A to a sparse matrix
970  MatConvert(A,
971  MATAIJ,
972 #if PETSC_VERSION_LESS_THAN(3, 7, 0)
973  MAT_REUSE_MATRIX,
974 #else
975  MAT_INPLACE_MATRIX,
976 #endif
977  &A);
978 
979  ISColoring iscoloring;
980  MatColoring mc;
981  MatColoringCreate(A, &mc);
982  MatColoringSetType(mc, coloring_algorithm);
983  MatColoringSetMaxColors(mc, static_cast<PetscInt>(colors));
984 
985  // Petsc normally colors by distance two (neighbors of neighbors), we just want one
986  MatColoringSetDistance(mc, 1);
987  MatColoringSetFromOptions(mc);
988  MatColoringApply(mc, &iscoloring);
989 
990  PetscInt nn;
991  IS * is;
992 #if PETSC_RELEASE_LESS_THAN(3, 12, 0)
993  ISColoringGetIS(iscoloring, &nn, &is);
994 #else
995  ISColoringGetIS(iscoloring, PETSC_USE_POINTER, &nn, &is);
996 #endif
997 
998  if (nn > static_cast<PetscInt>(colors))
999  throw std::runtime_error("Not able to color with designated number of colors");
1000 
1001  for (int i = 0; i < nn; i++)
1002  {
1003  PetscInt isize;
1004  const PetscInt * indices;
1005  ISGetLocalSize(is[i], &isize);
1006  ISGetIndices(is[i], &indices);
1007  for (int j = 0; j < isize; j++)
1008  {
1009  mooseAssert(indices[j] < static_cast<PetscInt>(vertex_colors.size()), "Index out of bounds");
1010  vertex_colors[indices[j]] = i;
1011  }
1012  ISRestoreIndices(is[i], &indices);
1013  }
1014 
1015  MatDestroy(&A);
1016  MatColoringDestroy(&mc);
1017  ISColoringDestroy(&iscoloring);
1018 }
PetscErrorCode PetscInt const PetscInt IS * is

◆ disableLinearConvergedReason()

void Moose::PetscSupport::disableLinearConvergedReason ( FEProblemBase fe_problem)

disable printing of the linear convergence reason

Definition at line 1034 of file PetscSupport.C.

Referenced by CommonOutputAction::act().

1035 {
1036  auto & petsc_options = fe_problem.getPetscOptions();
1037 
1038  petsc_options.flags.erase("-ksp_converged_reason");
1039 
1040  auto & pairs = petsc_options.pairs;
1041  auto it = MooseUtils::findPair(pairs, "-ksp_converged_reason", MooseUtils::Any);
1042  if (it != pairs.end())
1043  pairs.erase(it);
1044 }
Moose::PetscSupport::PetscOptions & getPetscOptions()
Retrieve a writable reference the PETSc options (used by PetscSupport)
MultiMooseEnum flags
Single value PETSc options (flags)
Definition: PetscSupport.h:47
C::iterator findPair(C &container, const M1 &first, const M2 &second)
Find a specific pair in a container matching on first, second or both pair components.
Definition: MooseUtils.h:1053
void erase(const std::string &names)
Un-assign a value.
static const struct MooseUtils::AnyType Any

◆ disableNonlinearConvergedReason()

void Moose::PetscSupport::disableNonlinearConvergedReason ( FEProblemBase fe_problem)

disable printing of the nonlinear convergence reason

Definition at line 1021 of file PetscSupport.C.

Referenced by CommonOutputAction::act().

1022 {
1023  auto & petsc_options = fe_problem.getPetscOptions();
1024 
1025  petsc_options.flags.erase("-snes_converged_reason");
1026 
1027  auto & pairs = petsc_options.pairs;
1028  auto it = MooseUtils::findPair(pairs, "-snes_converged_reason", MooseUtils::Any);
1029  if (it != pairs.end())
1030  pairs.erase(it);
1031 }
Moose::PetscSupport::PetscOptions & getPetscOptions()
Retrieve a writable reference the PETSc options (used by PetscSupport)
MultiMooseEnum flags
Single value PETSc options (flags)
Definition: PetscSupport.h:47
C::iterator findPair(C &container, const M1 &first, const M2 &second)
Find a specific pair in a container matching on first, second or both pair components.
Definition: MooseUtils.h:1053
void erase(const std::string &names)
Un-assign a value.
static const struct MooseUtils::AnyType Any

◆ getCommonPetscFlags()

MultiMooseEnum Moose::PetscSupport::getCommonPetscFlags ( )

A helper function to produce a MultiMooseEnum with commonly used PETSc single options (flags)

Definition at line 883 of file PetscSupport.C.

Referenced by getPetscValidParams(), AddFieldSplitAction::validParams(), and Split::validParams().

884 {
885  return MultiMooseEnum(
886  "-dm_moose_print_embedding -dm_view -ksp_converged_reason -ksp_gmres_modifiedgramschmidt "
887  "-ksp_monitor -ksp_monitor_snes_lg-snes_ksp_ew -ksp_snes_ew -snes_converged_reason "
888  "-snes_ksp -snes_ksp_ew -snes_linesearch_monitor -snes_mf -snes_mf_operator -snes_monitor "
889  "-snes_test_display -snes_view",
890  "",
891  true);
892 }
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...

◆ getCommonPetscKeys()

MultiMooseEnum Moose::PetscSupport::getCommonPetscKeys ( )

A helper function to produce a MultiMooseEnum with commonly used PETSc iname options (keys in key-value pairs)

Definition at line 895 of file PetscSupport.C.

Referenced by getPetscValidParams(), and AddFieldSplitAction::validParams().

896 {
897  return MultiMooseEnum("-ksp_atol -ksp_gmres_restart -ksp_max_it -ksp_pc_side -ksp_rtol "
898  "-ksp_type -mat_fd_coloring_err -mat_fd_type -mat_mffd_type "
899  "-pc_asm_overlap -pc_factor_levels "
900  "-pc_factor_mat_ordering_type -pc_hypre_boomeramg_grid_sweeps_all "
901  "-pc_hypre_boomeramg_max_iter "
902  "-pc_hypre_boomeramg_strong_threshold -pc_hypre_type -pc_type -snes_atol "
903  "-snes_linesearch_type "
904  "-snes_ls -snes_max_it -snes_rtol -snes_divergence_tolerance -snes_type "
905  "-sub_ksp_type -sub_pc_type",
906  "",
907  true);
908 }
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...

◆ getPetscKSPNormType()

KSPNormType Moose::PetscSupport::getPetscKSPNormType ( Moose::MooseKSPNormType  kspnorm)

Definition at line 433 of file PetscSupport.C.

Referenced by petscSetDefaultKSPNormType().

434 {
435  switch (kspnorm)
436  {
437  case Moose::KSPN_NONE:
438  return KSP_NORM_NONE;
440  return KSP_NORM_PRECONDITIONED;
442  return KSP_NORM_UNPRECONDITIONED;
443  case Moose::KSPN_NATURAL:
444  return KSP_NORM_NATURAL;
445  case Moose::KSPN_DEFAULT:
446  return KSP_NORM_DEFAULT;
447  default:
448  mooseError("Unknown KSP norm type requested.");
449  break;
450  }
451 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
Use whatever we have in PETSc.
Definition: MooseTypes.h:749

◆ getPetscPCSide()

PCSide Moose::PetscSupport::getPetscPCSide ( Moose::PCSideType  pcs)

Definition at line 416 of file PetscSupport.C.

Referenced by petscSetDefaultPCSide().

417 {
418  switch (pcs)
419  {
420  case Moose::PCS_LEFT:
421  return PC_LEFT;
422  case Moose::PCS_RIGHT:
423  return PC_RIGHT;
425  return PC_SYMMETRIC;
426  default:
427  mooseError("Unknown PC side requested.");
428  break;
429  }
430 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299

◆ getPetscValidLineSearches()

std::set< std::string > Moose::PetscSupport::getPetscValidLineSearches ( )

Returns the valid petsc line search options as a set of strings.

Definition at line 842 of file PetscSupport.C.

Referenced by FEProblemSolve::validParams().

843 {
844  return {"default", "shell", "none", "basic", "l2", "bt", "cp"};
845 }

◆ getPetscValidParams()

InputParameters Moose::PetscSupport::getPetscValidParams ( )

Returns the PETSc options that are common between Executioners and Preconditioners.

Returns
InputParameters object containing the PETSc related parameters

The output of this function should be added to the the parameters object of the overarching class

See also
CreateExecutionerAction

Definition at line 848 of file PetscSupport.C.

Referenced by FEProblemSolve::validParams(), and MoosePreconditioner::validParams().

849 {
851 
852  MooseEnum solve_type("PJFNK JFNK NEWTON FD LINEAR");
853  params.addParam<MooseEnum>("solve_type",
854  solve_type,
855  "PJFNK: Preconditioned Jacobian-Free Newton Krylov "
856  "JFNK: Jacobian-Free Newton Krylov "
857  "NEWTON: Full Newton Solve "
858  "FD: Use finite differences to compute Jacobian "
859  "LINEAR: Solving a linear problem");
860 
861  MooseEnum mffd_type("wp ds", "wp");
862  params.addParam<MooseEnum>("mffd_type",
863  mffd_type,
864  "Specifies the finite differencing type for "
865  "Jacobian-free solve types. Note that the "
866  "default is wp (for Walker and Pernice).");
867 
868  params.addParam<MultiMooseEnum>(
869  "petsc_options", getCommonPetscFlags(), "Singleton PETSc options");
870  params.addParam<MultiMooseEnum>(
871  "petsc_options_iname", getCommonPetscKeys(), "Names of PETSc name/value pairs");
872  params.addParam<std::vector<std::string>>(
873  "petsc_options_value",
874  "Values of PETSc name/value pairs (must correspond with \"petsc_options_iname\"");
875  params.addParamNamesToGroup("solve_type petsc_options petsc_options_iname petsc_options_value "
876  "mffd_type",
877  "PETSc");
878 
879  return params;
880 }
MultiMooseEnum getCommonPetscKeys()
A helper function to produce a MultiMooseEnum with commonly used PETSc iname options (keys in key-val...
Definition: PetscSupport.C:895
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
InputParameters emptyInputParameters()
MultiMooseEnum getCommonPetscFlags()
A helper function to produce a MultiMooseEnum with commonly used PETSc single options (flags) ...
Definition: PetscSupport.C:883
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:31
void addParam(const std::string &name, const S &value, const std::string &doc_string)
These methods add an option parameter and a documentation string to the InputParameters object...
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
void addParamNamesToGroup(const std::string &space_delim_names, const std::string group_name)
This method takes a space delimited list of parameter names and adds them to the specified group name...

◆ isSNESVI()

bool Moose::PetscSupport::isSNESVI ( FEProblemBase fe_problem)

check if SNES type is variational inequalities (VI) solver

Definition at line 911 of file PetscSupport.C.

Referenced by BoundsBase::BoundsBase().

912 {
913  PetscOptions & petsc = fe_problem.getPetscOptions();
914 
915  int argc;
916  char ** args;
917  PetscGetArgs(&argc, &args);
918 
919  std::vector<std::string> cml_arg;
920  for (int i = 0; i < argc; i++)
921  cml_arg.push_back(args[i]);
922 
923  if (MooseUtils::findPair(petsc.pairs, MooseUtils::Any, "vinewtonssls") == petsc.pairs.end() &&
924  MooseUtils::findPair(petsc.pairs, MooseUtils::Any, "vinewtonrsls") == petsc.pairs.end() &&
925  std::find(cml_arg.begin(), cml_arg.end(), "vinewtonssls") == cml_arg.end() &&
926  std::find(cml_arg.begin(), cml_arg.end(), "vinewtonrsls") == cml_arg.end())
927  return false;
928 
929  return true;
930 }
Moose::PetscSupport::PetscOptions & getPetscOptions()
Retrieve a writable reference the PETSc options (used by PetscSupport)
C::iterator findPair(C &container, const M1 &first, const M2 &second)
Find a specific pair in a container matching on first, second or both pair components.
Definition: MooseUtils.h:1053
static const struct MooseUtils::AnyType Any

◆ outputNorm()

void Moose::PetscSupport::outputNorm ( libMesh::Real  old_norm,
libMesh::Real  norm,
bool  use_color = false 
)

Helper function for outputting the norm values with/without color.

◆ petscLinearMonitor()

PetscErrorCode Moose::PetscSupport::petscLinearMonitor ( KSP  ,
PetscInt  its,
PetscReal  rnorm,
void void_ptr 
)

Helper function for displaying the linear residual during PETSC solve.

◆ petscNonlinearConverged()

PetscErrorCode Moose::PetscSupport::petscNonlinearConverged ( SNES  snes,
PetscInt  it,
PetscReal  xnorm,
PetscReal  snorm,
PetscReal  fnorm,
SNESConvergedReason *  reason,
void ctx 
)

Definition at line 277 of file PetscSupport.C.

Referenced by petscSetDefaults().

284 {
285  FEProblemBase & problem = *static_cast<FEProblemBase *>(ctx);
286 
287  // Let's be nice and always check PETSc error codes.
288  PetscErrorCode ierr = 0;
289 
290  // Temporary variables to store SNES tolerances. Usual C-style would be to declare
291  // but not initialize these... but it bothers me to leave anything uninitialized.
292  PetscReal atol = 0.; // absolute convergence tolerance
293  PetscReal rtol = 0.; // relative convergence tolerance
294  PetscReal stol = 0.; // convergence (step) tolerance in terms of the norm of the change in the
295  // solution between steps
296  PetscInt maxit = 0; // maximum number of iterations
297  PetscInt maxf = 0; // maximum number of function evaluations
298 
299  // Ask the SNES object about its tolerances.
300  ierr = SNESGetTolerances(snes, &atol, &rtol, &stol, &maxit, &maxf);
301  CHKERRABORT(problem.comm().get(), ierr);
302 
303  // Ask the SNES object about its divergence tolerance.
304  PetscReal divtol = 0.; // relative divergence tolerance
305 #if !PETSC_VERSION_LESS_THAN(3, 8, 0)
306  ierr = SNESGetDivergenceTolerance(snes, &divtol);
307  CHKERRABORT(problem.comm().get(), ierr);
308 #endif
309 
310  // Get current number of function evaluations done by SNES.
311  PetscInt nfuncs = 0;
312  ierr = SNESGetNumberFunctionEvals(snes, &nfuncs);
313  CHKERRABORT(problem.comm().get(), ierr);
314 
315  // Whether or not to force SNESSolve() take at least one iteration regardless of the initial
316  // residual norm
317 #if !PETSC_VERSION_LESS_THAN(3, 8, 4)
318  PetscBool force_iteration = PETSC_FALSE;
319  ierr = SNESGetForceIteration(snes, &force_iteration);
320  CHKERRABORT(problem.comm().get(), ierr);
321 
322  if (force_iteration && !(problem.getNonlinearForcedIterations()))
323  problem.setNonlinearForcedIterations(1);
324 
325  if (!force_iteration && (problem.getNonlinearForcedIterations()))
326  {
327  ierr = SNESSetForceIteration(snes, PETSC_TRUE);
328  CHKERRABORT(problem.comm().get(), ierr);
329  }
330 #endif
331 
332  // See if SNESSetFunctionDomainError() has been called. Note:
333  // SNESSetFunctionDomainError() and SNESGetFunctionDomainError()
334  // were added in different releases of PETSc.
335  PetscBool domainerror;
336  ierr = SNESGetFunctionDomainError(snes, &domainerror);
337  CHKERRABORT(problem.comm().get(), ierr);
338  if (domainerror)
339  {
340  *reason = SNES_DIVERGED_FUNCTION_DOMAIN;
341  return 0;
342  }
343 
344  // Error message that will be set by the FEProblemBase.
345  std::string msg;
346 
347  // xnorm: 2-norm of current iterate
348  // snorm: 2-norm of current step
349  // fnorm: 2-norm of function at current iterate
350  MooseNonlinearConvergenceReason moose_reason =
351  problem.checkNonlinearConvergence(msg,
352  it,
353  xnorm,
354  snorm,
355  fnorm,
356  rtol,
357  divtol,
358  stol,
359  atol,
360  nfuncs,
361  maxf,
363 
364  if (msg.length() > 0)
365 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
366  PetscInfo(snes, "%s", msg.c_str());
367 #else
368  PetscInfo(snes, msg.c_str());
369 #endif
370 
371  switch (moose_reason)
372  {
374  *reason = SNES_CONVERGED_ITERATING;
375  break;
376 
378  *reason = SNES_CONVERGED_FNORM_ABS;
379  break;
380 
382  *reason = SNES_CONVERGED_FNORM_RELATIVE;
383  break;
384 
386 #if !PETSC_VERSION_LESS_THAN(3, 8, 0) // A new convergence enum in PETSc 3.8
387  *reason = SNES_DIVERGED_DTOL;
388 #endif
389  break;
390 
392  *reason = SNES_CONVERGED_SNORM_RELATIVE;
393  break;
394 
396  *reason = SNES_DIVERGED_FUNCTION_COUNT;
397  break;
398 
400  *reason = SNES_DIVERGED_FNORM_NAN;
401  break;
402 
404  *reason = SNES_DIVERGED_LINE_SEARCH;
405  break;
406 
408  *reason = SNES_DIVERGED_LOCAL_MIN;
409  break;
410  }
411 
412  return 0;
413 }
const Parallel::Communicator & comm() const
virtual MooseNonlinearConvergenceReason checkNonlinearConvergence(std::string &msg, const PetscInt it, const Real xnorm, const Real snorm, const Real fnorm, const Real rtol, const Real divtol, const Real stol, const Real abstol, const PetscInt nfuncs, const PetscInt max_funcs, const Real div_threshold)
Check for convergence of the nonlinear solution.
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
ierr
auto max(const L &left, const R &right)
void setNonlinearForcedIterations(const unsigned int nl_forced_its)
method setting the minimum number of nonlinear iterations before performing divergence checks ...
MooseNonlinearConvergenceReason
Enumeration for nonlinear convergence reasons.
void * ctx
unsigned int getNonlinearForcedIterations() const
method returning the number of forced nonlinear iterations

◆ petscSetDefaultKSPNormType()

void Moose::PetscSupport::petscSetDefaultKSPNormType ( FEProblemBase problem,
KSP  ksp 
)

Set norm type.

Definition at line 454 of file PetscSupport.C.

Referenced by Moose::SlepcSupport::mooseSlepcEPSSNESKSPSetPCSide(), and petscSetKSPDefaults().

455 {
456  for (const auto i : make_range(problem.numSolverSystems()))
457  {
458  SolverSystem & sys = problem.getSolverSystem(i);
459  KSPSetNormType(ksp, getPetscKSPNormType(sys.getMooseKSPNormType()));
460  }
461 }
KSPNormType getPetscKSPNormType(Moose::MooseKSPNormType kspnorm)
Definition: PetscSupport.C:433
IntRange< T > make_range(T beg, T end)
SolverSystem & getSolverSystem(unsigned int sys_num)
Get non-constant reference to a solver system.
Moose::MooseKSPNormType getMooseKSPNormType()
Get the norm in which the linear convergence is measured.
Definition: SolverSystem.h:67
virtual std::size_t numSolverSystems() const override

◆ petscSetDefaultPCSide()

void Moose::PetscSupport::petscSetDefaultPCSide ( FEProblemBase problem,
KSP  ksp 
)

Setup which side we want to apply preconditioner.

Definition at line 464 of file PetscSupport.C.

Referenced by Moose::SlepcSupport::mooseSlepcEPSSNESKSPSetPCSide(), and petscSetKSPDefaults().

465 {
466  for (const auto i : make_range(problem.numSolverSystems()))
467  {
468  SolverSystem & sys = problem.getSolverSystem(i);
469 
470  // PETSc 3.2.x+
471  if (sys.getPCSide() != Moose::PCS_DEFAULT)
472  KSPSetPCSide(ksp, getPetscPCSide(sys.getPCSide()));
473  }
474 }
PCSide getPetscPCSide(Moose::PCSideType pcs)
Definition: PetscSupport.C:416
Use whatever we have in PETSc.
Definition: MooseTypes.h:737
IntRange< T > make_range(T beg, T end)
SolverSystem & getSolverSystem(unsigned int sys_num)
Get non-constant reference to a solver system.
virtual std::size_t numSolverSystems() const override
Moose::PCSideType getPCSide()
Get the current preconditioner side.
Definition: SolverSystem.h:56

◆ petscSetDefaults()

void Moose::PetscSupport::petscSetDefaults ( FEProblemBase problem)

Sets the default options for PETSc.

Definition at line 499 of file PetscSupport.C.

Referenced by FEProblemBase::initPetscOutputAndSomeSolverSettings(), and Moose::setSolverDefaults().

500 {
501  for (auto nl_index : make_range(problem.numNonlinearSystems()))
502  {
503  // dig out PETSc solver
504  NonlinearSystemBase & nl = problem.getNonlinearSystemBase(nl_index);
505  PetscNonlinearSolver<Number> * petsc_solver =
506  dynamic_cast<PetscNonlinearSolver<Number> *>(nl.nonlinearSolver());
507  SNES snes = petsc_solver->snes();
508  KSP ksp;
509  auto ierr = SNESGetKSP(snes, &ksp);
510  CHKERRABORT(nl.comm().get(), ierr);
511 
512  ierr = SNESSetMaxLinearSolveFailures(snes, 1000000);
513  CHKERRABORT(nl.comm().get(), ierr);
514 
515  ierr = SNESSetCheckJacobianDomainError(snes, PETSC_TRUE);
516  CHKERRABORT(nl.comm().get(), ierr);
517 
518  // In 3.0.0, the context pointer must actually be used, and the
519  // final argument to KSPSetConvergenceTest() is a pointer to a
520  // routine for destroying said private data context. In this case,
521  // we use the default context provided by PETSc in addition to
522  // a few other tests.
523  {
524  ierr = SNESSetConvergenceTest(snes, petscNonlinearConverged, &problem, LIBMESH_PETSC_NULLPTR);
525  CHKERRABORT(nl.comm().get(), ierr);
526  }
527 
528  petscSetKSPDefaults(problem, ksp);
529  }
530 
531  for (auto sys_index : make_range(problem.numLinearSystems()))
532  {
533  // dig out PETSc solver
534  LinearSystem & lin_sys = problem.getLinearSystem(sys_index);
535  PetscLinearSolver<Number> * petsc_solver = dynamic_cast<PetscLinearSolver<Number> *>(
536  lin_sys.linearImplicitSystem().get_linear_solver());
537  KSP ksp = petsc_solver->ksp();
538  petscSetKSPDefaults(problem, ksp);
539  }
540 }
virtual NonlinearSolver< Number > * nonlinearSolver()=0
virtual std::size_t numNonlinearSystems() const override
const Parallel::Communicator & comm() const
void petscSetKSPDefaults(FEProblemBase &problem, KSP ksp)
Set the default options for a KSP.
Definition: PetscSupport.C:477
ierr
Nonlinear system to be solved.
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
LinearSystem & getLinearSystem(unsigned int sys_num)
Get non-constant reference to a linear system.
IntRange< T > make_range(T beg, T end)
Linear system to be solved.
Definition: LinearSystem.h:35
virtual std::size_t numLinearSystems() const override
PetscErrorCode petscNonlinearConverged(SNES snes, PetscInt it, PetscReal xnorm, PetscReal snorm, PetscReal fnorm, SNESConvergedReason *reason, void *ctx)
Definition: PetscSupport.C:277
LinearImplicitSystem & linearImplicitSystem()
Return a reference to the stored linear implicit system.
Definition: LinearSystem.h:72

◆ petscSetKSPDefaults()

void Moose::PetscSupport::petscSetKSPDefaults ( FEProblemBase problem,
KSP  ksp 
)

Set the default options for a KSP.

Definition at line 477 of file PetscSupport.C.

Referenced by petscSetDefaults(), and setLinearSolverDefaults().

478 {
479  auto & es = problem.es();
480 
481  PetscReal rtol = es.parameters.get<Real>("linear solver tolerance");
482  PetscReal atol = es.parameters.get<Real>("linear solver absolute tolerance");
483 
484  // MOOSE defaults this to -1 for some dumb reason
485  if (atol < 0)
486  atol = 1e-50;
487 
488  PetscReal maxits = es.parameters.get<unsigned int>("linear solver maximum iterations");
489 
490  // 1e100 is because we don't use divtol currently
491  KSPSetTolerances(ksp, rtol, atol, 1e100, maxits);
492 
493  petscSetDefaultPCSide(problem, ksp);
494 
495  petscSetDefaultKSPNormType(problem, ksp);
496 }
void petscSetDefaultKSPNormType(FEProblemBase &problem, KSP ksp)
Set norm type.
Definition: PetscSupport.C:454
virtual EquationSystems & es() override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
void petscSetDefaultPCSide(FEProblemBase &problem, KSP ksp)
Setup which side we want to apply preconditioner.
Definition: PetscSupport.C:464

◆ petscSetOptions()

void Moose::PetscSupport::petscSetOptions ( const PetscOptions po,
const SolverParams solver_params 
)

A function for setting the PETSc options in PETSc from the options supplied to MOOSE.

Definition at line 239 of file PetscSupport.C.

Referenced by Moose::SlepcSupport::slepcSetOptions(), FEProblemBase::solve(), and FEProblemBase::solveLinearSystem().

240 {
241 #if PETSC_VERSION_LESS_THAN(3, 7, 0)
242  PetscOptionsClear();
243 #else
244  PetscOptionsClear(LIBMESH_PETSC_NULLPTR);
245 #endif
246 
247  setSolverOptions(solver_params);
248 
249  // Add any additional options specified in the input file
250  for (const auto & flag : po.flags)
251  setSinglePetscOption(flag.rawName().c_str());
252 
253  // Add option pairs
254  for (auto & option : po.pairs)
255  setSinglePetscOption(option.first, option.second);
256 
258 }
void addPetscOptionsFromCommandline()
Definition: PetscSupport.C:221
void setSolverOptions(const SolverParams &solver_params)
Definition: PetscSupport.C:139
void setSinglePetscOption(const std::string &name, const std::string &value="")
A wrapper function for dealing with different versions of PetscOptionsSetValue.
Definition: PetscSupport.C:933

◆ petscSetupDM()

void Moose::PetscSupport::petscSetupDM ( NonlinearSystemBase nl,
const std::string &  dm_name 
)

Setup the PETSc DM object.

Definition at line 176 of file PetscSupport.C.

Referenced by NonlinearSystemBase::setupDM().

177 {
178  PetscErrorCode ierr;
179  PetscBool ismoose;
180  DM dm = LIBMESH_PETSC_NULLPTR;
181 
182  // Initialize the part of the DM package that's packaged with Moose; in the PETSc source tree this
183  // call would be in DMInitializePackage()
185  CHKERRABORT(nl.comm().get(), ierr);
186  // Create and set up the DM that will consume the split options and deal with block matrices.
187  PetscNonlinearSolver<Number> * petsc_solver =
188  dynamic_cast<PetscNonlinearSolver<Number> *>(nl.nonlinearSolver());
189  SNES snes = petsc_solver->snes();
190  // if there exists a DMMoose object, not to recreate a new one
191  ierr = SNESGetDM(snes, &dm);
192  CHKERRABORT(nl.comm().get(), ierr);
193  if (dm)
194  {
195  ierr = PetscObjectTypeCompare((PetscObject)dm, DMMOOSE, &ismoose);
196  CHKERRABORT(nl.comm().get(), ierr);
197  if (ismoose)
198  return;
199  }
200  ierr = DMCreateMoose(nl.comm().get(), nl, dm_name, &dm);
201  CHKERRABORT(nl.comm().get(), ierr);
202  ierr = DMSetFromOptions(dm);
203  CHKERRABORT(nl.comm().get(), ierr);
204  ierr = DMSetUp(dm);
205  CHKERRABORT(nl.comm().get(), ierr);
206  ierr = SNESSetDM(snes, dm);
207  CHKERRABORT(nl.comm().get(), ierr);
208  ierr = DMDestroy(&dm);
209  CHKERRABORT(nl.comm().get(), ierr);
210  // We temporarily comment out this updating function because
211  // we lack an approach to check if the problem
212  // structure has been changed from the last iteration.
213  // The indices will be rebuilt for every timestep.
214  // TODO: figure out a way to check the structure changes of the
215  // matrix
216  // ierr = SNESSetUpdate(snes,SNESUpdateDMMoose);
217  // CHKERRABORT(nl.comm().get(),ierr);
218 }
virtual NonlinearSolver< Number > * nonlinearSolver()=0
const Parallel::Communicator & comm() const
ierr
PetscErrorCode DMCreateMoose(MPI_Comm comm, NonlinearSystemBase &nl, const std::string &dm_name, DM *dm)
Create a MOOSE DM.
PetscErrorCode DMMooseRegisterAll()

◆ petscSetupOutput()

PetscErrorCode Moose::PetscSupport::petscSetupOutput ( CommandLine cmd_line)

Definition at line 261 of file PetscSupport.C.

Referenced by MooseApp::executeExecutioner().

262 {
263  char code[10] = {45, 45, 109, 111, 111, 115, 101};
264  const std::vector<std::string> argv = cmd_line->getArguments();
265  for (const auto & arg : argv)
266  {
267  if (arg.compare(code) == 0)
268  {
270  break;
271  }
272  }
273  return 0;
274 }
static void petscSetupOutput()
Output string for setting up PETSC output.
Definition: Console.C:810
const std::vector< std::string > & getArguments()
Return the raw argv arguments as a vector.
Definition: CommandLine.h:71

◆ processPetscFlags()

void Moose::PetscSupport::processPetscFlags ( const MultiMooseEnum petsc_flags,
PetscOptions petsc_options 
)

Populate flags in a given PetscOptions object using a vector of input arguments.

Parameters
petsc_flagsContainer holding the flags of the petsc options
petsc_optionsData structure which handles petsc options within moose

"-log_summary" cannot be used in the input file. This option needs to be set when PETSc is initialized which happens before the parser is even created. We'll throw an error if somebody attempts to add this option later.

"-log_summary" cannot be used in the input file. This option needs to be set when PETSc is initialized which happens before the parser is even created. We'll throw an error if somebody attempts to add this option later.

Definition at line 603 of file PetscSupport.C.

Referenced by storePetscOptions().

604 {
605  // Update the PETSc single flags
606  for (const auto & option : petsc_flags)
607  {
614  if (option == "-log_summary" || option == "-log_view")
615  mooseError("The PETSc option \"-log_summary\" or \"-log_view\" can only be used on the "
616  "command line. Please "
617  "remove it from the input file");
618 
619  // Warn about superseded PETSc options (Note: -snes is not a REAL option, but people used it in
620  // their input files)
621  else
622  {
623  std::string help_string;
624  if (option == "-snes" || option == "-snes_mf" || option == "-snes_mf_operator")
625  help_string = "Please set the solver type through \"solve_type\".";
626  else if (option == "-ksp_monitor")
627  help_string = "Please use \"Outputs/print_linear_residuals=true\"";
628 
629  if (help_string != "")
630  mooseWarning("The PETSc option ",
631  std::string(option),
632  " should not be used directly in a MOOSE input file. ",
633  help_string);
634  }
635 
636  // Update the stored items, but do not create duplicates
637  if (!po.flags.contains(option))
638  po.flags.push_back(option);
639  }
640 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
void mooseWarning(Args &&... args)
Emit a warning message with the given stringified, concatenated args.
Definition: MooseError.h:333

◆ processPetscPairs()

void Moose::PetscSupport::processPetscPairs ( const std::vector< std::pair< MooseEnumItem, std::string >> &  petsc_pair_options,
const unsigned int  mesh_dimension,
PetscOptions petsc_options 
)

Populate name and value pairs in a given PetscOptions object using vectors of input arguments.

Parameters
petsc_pair_optionsOption-value pairs of petsc settings
mesh_dimensionThe mesh dimension, needed for multigrid settings
petsc_optionsData structure which handles petsc options within moose

Definition at line 643 of file PetscSupport.C.

Referenced by DiffusionPhysicsBase::addPreconditioning(), and storePetscOptions().

646 {
647  // the boolean in these pairs denote whether the user has specified any of the reason flags in the
648  // input file
649  std::array<std::pair<bool, std::string>, 2> reason_flags = {
650  {std::make_pair(false, "-snes_converged_reason"),
651  std::make_pair(false, "-ksp_converged_reason")}};
652 
653  for (auto & reason_flag : reason_flags)
654  if (po.flags.contains(reason_flag.second))
655  // We register the reason option as already existing
656  reason_flag.first = true;
657 
658  // Setup the name value pairs
659  bool boomeramg_found = false;
660  bool strong_threshold_found = false;
661 #if !PETSC_VERSION_LESS_THAN(3, 7, 0)
662  bool superlu_dist_found = false;
663  bool fact_pattern_found = false;
664  bool tiny_pivot_found = false;
665 #endif
666  std::string pc_description = "";
667 #if !PETSC_VERSION_LESS_THAN(3, 12, 0)
668  // If users use HMG, we would like to set
669  bool hmg_found = false;
670  bool matptap_found = false;
671  bool hmg_strong_threshold_found = false;
672 #endif
673 #if !PETSC_VERSION_LESS_THAN(3, 19, 2)
674  // Check for if users have set the options_left option
675  bool options_left_set = false;
676 #endif
677  std::vector<std::pair<std::string, std::string>> new_options;
678 
679  for (const auto & option : petsc_pair_options)
680  {
681  new_options.clear();
682 
683  // Do not add duplicate settings
684  if (MooseUtils::findPair(po.pairs, option.first, MooseUtils::Any) == po.pairs.end())
685  {
686 #if !PETSC_VERSION_LESS_THAN(3, 9, 0)
687  if (option.first == "-pc_factor_mat_solver_package")
688  new_options.emplace_back("-pc_factor_mat_solver_type", option.second);
689 #else
690  if (option.first == "-pc_factor_mat_solver_type")
691  new_options.push_back("-pc_factor_mat_solver_package", option.second);
692 #endif
693 
694  // Look for a pc description
695  if (option.first == "-pc_type" || option.first == "-pc_sub_type" ||
696  option.first == "-pc_hypre_type")
697  pc_description += option.second + ' ';
698 
699 #if !PETSC_VERSION_LESS_THAN(3, 12, 0)
700  if (option.first == "-pc_type" && option.second == "hmg")
701  hmg_found = true;
702 
703  // MPIAIJ for PETSc 3.12.0: -matptap_via
704  // MAIJ for PETSc 3.12.0: -matmaijptap_via
705  // MPIAIJ for PETSc 3.13 to 3.16: -matptap_via, -matproduct_ptap_via
706  // MAIJ for PETSc 3.13 to 3.16: -matproduct_ptap_via
707  // MPIAIJ for PETSc 3.17 and higher: -matptap_via, -mat_product_algorithm
708  // MAIJ for PETSc 3.17 and higher: -mat_product_algorithm
709 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
710  if (hmg_found && (option.first == "-matptap_via" || option.first == "-matmaijptap_via" ||
711  option.first == "-matproduct_ptap_via"))
712  new_options.emplace_back("-mat_product_algorithm", option.second);
713 #elif !PETSC_VERSION_LESS_THAN(3, 13, 0)
714  if (hmg_found && (option.first == "-matptap_via" || option.first == "-matmaijptap_via"))
715  new_options.emplace_back("-matproduct_ptap_via", option.second);
716 #else
717  if (hmg_found && (option.first == "-matproduct_ptap_via"))
718  {
719  new_options.emplace_back("-matptap_via", option.second);
720  new_options.emplace_back("-matmaijptap_via", option.second);
721  }
722 #endif
723 
724  if (option.first == "-matptap_via" || option.first == "-matmaijptap_via" ||
725  option.first == "-matproduct_ptap_via" || option.first == "-mat_product_algorithm")
726  matptap_found = true;
727 
728  // For 3D problems, we need to set this 0.7
729  if (option.first == "-hmg_inner_pc_hypre_boomeramg_strong_threshold")
730  hmg_strong_threshold_found = true;
731 #endif
732  // This special case is common enough that we'd like to handle it for the user.
733  if (option.first == "-pc_hypre_type" && option.second == "boomeramg")
734  boomeramg_found = true;
735  if (option.first == "-pc_hypre_boomeramg_strong_threshold")
736  strong_threshold_found = true;
737 #if !PETSC_VERSION_LESS_THAN(3, 7, 0)
738  if ((option.first == "-pc_factor_mat_solver_package" ||
739  option.first == "-pc_factor_mat_solver_type") &&
740  option.second == "superlu_dist")
741  superlu_dist_found = true;
742  if (option.first == "-mat_superlu_dist_fact")
743  fact_pattern_found = true;
744  if (option.first == "-mat_superlu_dist_replacetinypivot")
745  tiny_pivot_found = true;
746 #endif
747 
748 #if !PETSC_VERSION_LESS_THAN(3, 19, 2)
749  if (option.first == "-options_left")
750  options_left_set = true;
751 #endif
752 
753  if (!new_options.empty())
754  std::copy(new_options.begin(), new_options.end(), std::back_inserter(po.pairs));
755  else
756  po.pairs.push_back(option);
757  }
758  else
759  {
760  for (unsigned int j = 0; j < po.pairs.size(); j++)
761  if (option.first == po.pairs[j].first)
762  po.pairs[j].second = option.second;
763  }
764  }
765 
766 #if !PETSC_VERSION_LESS_THAN(3, 14, 0)
767  for (const auto & reason_flag : reason_flags)
768  // Was the option already found in PetscOptions::flags? Or does it exist in PetscOptions::pairs
769  // as an iname already? If not, then we add our flag
770  if (!reason_flag.first && (std::find_if(po.pairs.begin(),
771  po.pairs.end(),
772  [&reason_flag](auto & pair) {
773  return pair.first == reason_flag.second;
774  }) == po.pairs.end()))
775  po.pairs.emplace_back(reason_flag.second, "::failed");
776 #endif
777 
778  // When running a 3D mesh with boomeramg, it is almost always best to supply a strong threshold
779  // value. We will provide that for the user here if they haven't supplied it themselves.
780  if (boomeramg_found && !strong_threshold_found && mesh_dimension == 3)
781  {
782  po.pairs.emplace_back("-pc_hypre_boomeramg_strong_threshold", "0.7");
783  pc_description += "strong_threshold: 0.7 (auto)";
784  }
785 
786 #if !PETSC_VERSION_LESS_THAN(3, 12, 0)
787  if (hmg_found && !hmg_strong_threshold_found && mesh_dimension == 3)
788  {
789  po.pairs.emplace_back("-hmg_inner_pc_hypre_boomeramg_strong_threshold", "0.7");
790  pc_description += "strong_threshold: 0.7 (auto)";
791  }
792 
793  // Default PETSc PtAP takes too much memory, and it is not quite useful
794  // Let us switch to use new algorithm
795  if (hmg_found && !matptap_found)
796  {
797 #if !PETSC_VERSION_LESS_THAN(3, 17, 0)
798  po.pairs.emplace_back("-mat_product_algorithm", "allatonce");
799 #elif !PETSC_VERSION_LESS_THAN(3, 13, 0)
800  po.pairs.emplace_back("-matproduct_ptap_via", "allatonce");
801 #else
802  po.pairs.emplace_back("-matptap_via", "allatonce");
803  po.pairs.emplace_back("-matmaijptap_via", "allatonce");
804 #endif
805  }
806 #endif
807 
808 #if !PETSC_VERSION_LESS_THAN(3, 7, 0)
809  // In PETSc-3.7.{0--4}, there is a bug when using superlu_dist, and we have to use
810  // SamePattern_SameRowPerm, otherwise we use whatever we have in PETSc
811  if (superlu_dist_found && !fact_pattern_found)
812  {
813  po.pairs.emplace_back("-mat_superlu_dist_fact",
814 #if PETSC_VERSION_LESS_THAN(3, 7, 5)
815  "SamePattern_SameRowPerm");
816  pc_description += "mat_superlu_dist_fact: SamePattern_SameRowPerm ";
817 #else
818  "SamePattern");
819  pc_description += "mat_superlu_dist_fact: SamePattern ";
820 #endif
821  }
822 
823  // restore this superlu option
824  if (superlu_dist_found && !tiny_pivot_found)
825  {
826  po.pairs.emplace_back("-mat_superlu_dist_replacetinypivot", "1");
827  pc_description += " mat_superlu_dist_replacetinypivot: true ";
828  }
829 #endif
830  // Set Preconditioner description
831  po.pc_description += pc_description;
832 
833  // Turn off default options_left warnings added in 3.19.3 pre-release for all PETSc builds
834  // (PETSc commit: 59f199a7), unless the user has set a preference.
835 #if !PETSC_VERSION_LESS_THAN(3, 19, 2)
836  if (!options_left_set && !po.flags.contains("-options_left"))
837  po.pairs.emplace_back("-options_left", "0");
838 #endif
839 }
C::iterator findPair(C &container, const M1 &first, const M2 &second)
Find a specific pair in a container matching on first, second or both pair components.
Definition: MooseUtils.h:1053
static const struct MooseUtils::AnyType Any

◆ setLinearSolverDefaults()

template<typename T >
void Moose::PetscSupport::setLinearSolverDefaults ( FEProblemBase problem,
LinearSolver< T > &  linear_solver 
)

Set the defaults for a libMesh LinearSolver.

Used in explicit solves

Definition at line 70 of file PetscSupport.h.

Referenced by ExplicitTimeIntegrator::meshChanged().

71 {
72  petscSetKSPDefaults(problem, libMesh::cast_ref<PetscLinearSolver<T> &>(linear_solver).ksp());
73 }
Tnew cast_ref(Told &oldvar)
void petscSetKSPDefaults(FEProblemBase &problem, KSP ksp)
Set the default options for a KSP.
Definition: PetscSupport.C:477

◆ setSinglePetscOption()

void Moose::PetscSupport::setSinglePetscOption ( const std::string &  name,
const std::string &  value = "" 
)

A wrapper function for dealing with different versions of PetscOptionsSetValue.

This is not generally called from MOOSE code, it is instead intended to be called by stuff in MOOSE::PetscSupport.

Definition at line 933 of file PetscSupport.C.

Referenced by Moose::SlepcSupport::clearFreeNonlinearPowerIterations(), petscSetOptions(), Moose::SlepcSupport::setEigenProblemOptions(), Moose::SlepcSupport::setEigenSolverOptions(), Moose::SlepcSupport::setFreeNonlinearPowerIterations(), Moose::SlepcSupport::setNewtonPetscOptions(), Moose::SlepcSupport::setNonlinearPowerOptions(), Moose::SlepcSupport::setSlepcEigenSolverTolerances(), setSolverOptions(), Moose::SlepcSupport::setWhichEigenPairsOptions(), and Moose::SlepcSupport::slepcSetOptions().

934 {
935  PetscErrorCode ierr;
936 
937 #if PETSC_VERSION_LESS_THAN(3, 7, 0)
938  ierr = PetscOptionsSetValue(name.c_str(), value == "" ? LIBMESH_PETSC_NULLPTR : value.c_str());
939 #else
940  // PETSc 3.7.0 and later version. First argument is the options
941  // database to use, NULL indicates the default global database.
942  ierr = PetscOptionsSetValue(
943  LIBMESH_PETSC_NULLPTR, name.c_str(), value == "" ? LIBMESH_PETSC_NULLPTR : value.c_str());
944 #endif
945 
946  // Not convenient to use the usual error checking macro, because we
947  // don't have a specific communicator in this helper function.
948  if (ierr)
949  mooseError("Error setting PETSc option: ", name);
950 }
std::string name(const ElemQuality q)
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
ierr
Real value(unsigned n, unsigned alpha, unsigned beta, Real x)

◆ setSolverOptions()

void Moose::PetscSupport::setSolverOptions ( const SolverParams solver_params)

Definition at line 139 of file PetscSupport.C.

Referenced by petscSetOptions().

140 {
141  // set PETSc options implied by a solve type
142  switch (solver_params._type)
143  {
144  case Moose::ST_PJFNK:
145  setSinglePetscOption("-snes_mf_operator");
146  setSinglePetscOption("-mat_mffd_type", stringify(solver_params._mffd_type));
147  break;
148 
149  case Moose::ST_JFNK:
150  setSinglePetscOption("-snes_mf");
151  setSinglePetscOption("-mat_mffd_type", stringify(solver_params._mffd_type));
152  break;
153 
154  case Moose::ST_NEWTON:
155  break;
156 
157  case Moose::ST_FD:
158  setSinglePetscOption("-snes_fd");
159  break;
160 
161  case Moose::ST_LINEAR:
162  setSinglePetscOption("-snes_type", "ksponly");
163  setSinglePetscOption("-snes_monitor_cancel");
164  break;
165  }
166 
167  Moose::LineSearchType ls_type = solver_params._line_search;
168  if (ls_type == Moose::LS_NONE)
169  ls_type = Moose::LS_BASIC;
170 
171  if (ls_type != Moose::LS_DEFAULT && ls_type != Moose::LS_CONTACT && ls_type != Moose::LS_PROJECT)
172  setSinglePetscOption("-snes_linesearch_type", stringify(ls_type));
173 }
Full Newton Solve.
Definition: MooseTypes.h:759
Moose::LineSearchType _line_search
Definition: SolverParams.h:20
Solving a linear problem.
Definition: MooseTypes.h:761
Moose::MffdType _mffd_type
Definition: SolverParams.h:21
LineSearchType
Type of the line search.
Definition: MooseTypes.h:838
std::string stringify(const MffdType &t)
Definition: PetscSupport.C:124
Jacobian-Free Newton Krylov.
Definition: MooseTypes.h:758
Moose::SolveType _type
Definition: SolverParams.h:19
Use finite differences to compute Jacobian.
Definition: MooseTypes.h:760
Preconditioned Jacobian-Free Newton Krylov.
Definition: MooseTypes.h:757
void setSinglePetscOption(const std::string &name, const std::string &value="")
A wrapper function for dealing with different versions of PetscOptionsSetValue.
Definition: PetscSupport.C:933

◆ storePetscOptions()

void Moose::PetscSupport::storePetscOptions ( FEProblemBase fe_problem,
const InputParameters params 
)

Stores the PETSc options supplied from the InputParameters with MOOSE.

Definition at line 543 of file PetscSupport.C.

Referenced by SetupPreconditionerAction::act(), and FEProblemSolve::FEProblemSolve().

544 {
545  // Note: Options set in the Preconditioner block will override those set in the Executioner block
546  if (params.isParamValid("solve_type") && !params.isParamValid("_use_eigen_value"))
547  {
548  // Extract the solve type
549  const std::string & solve_type = params.get<MooseEnum>("solve_type");
550  fe_problem.solverParams()._type = Moose::stringToEnum<Moose::SolveType>(solve_type);
551  }
552 
553  if (params.isParamValid("line_search"))
554  {
555  MooseEnum line_search = params.get<MooseEnum>("line_search");
556  if (fe_problem.solverParams()._line_search == Moose::LS_INVALID || line_search != "default")
557  {
558  Moose::LineSearchType enum_line_search =
559  Moose::stringToEnum<Moose::LineSearchType>(line_search);
560  fe_problem.solverParams()._line_search = enum_line_search;
561  if (enum_line_search == LS_CONTACT || enum_line_search == LS_PROJECT)
562  for (auto nl_index : make_range(fe_problem.numNonlinearSystems()))
563  {
564  NonlinearImplicitSystem * nl_system = dynamic_cast<NonlinearImplicitSystem *>(
565  &fe_problem.getNonlinearSystemBase(nl_index).system());
566  if (!nl_system)
567  mooseError("You've requested a line search but you must be solving an EigenProblem. "
568  "These two things are not consistent.");
569  PetscNonlinearSolver<Real> * petsc_nonlinear_solver =
570  dynamic_cast<PetscNonlinearSolver<Real> *>(nl_system->nonlinear_solver.get());
571  if (!petsc_nonlinear_solver)
572  mooseError("Currently the MOOSE line searches all use Petsc, so you "
573  "must use Petsc as your non-linear solver.");
574  petsc_nonlinear_solver->linesearch_object =
575  std::make_unique<ComputeLineSearchObjectWrapper>(fe_problem);
576  }
577  }
578  }
579 
580  if (params.isParamValid("mffd_type"))
581  {
582  MooseEnum mffd_type = params.get<MooseEnum>("mffd_type");
583  fe_problem.solverParams()._mffd_type = Moose::stringToEnum<Moose::MffdType>(mffd_type);
584  }
585 
586  // The parameters contained in the Action
587  const auto & petsc_options = params.get<MultiMooseEnum>("petsc_options");
588  const auto & petsc_pair_options =
589  params.get<MooseEnumItem, std::string>("petsc_options_iname", "petsc_options_value");
590 
591  // A reference to the PetscOptions object that contains the settings that will be used in the
592  // solve
594 
595  // First process the single petsc options/flags
596  processPetscFlags(petsc_options, po);
597 
598  // Then process the option-value pairs
599  processPetscPairs(petsc_pair_options, fe_problem.mesh().dimension(), po);
600 }
Moose::PetscSupport::PetscOptions & getPetscOptions()
Retrieve a writable reference the PETSc options (used by PetscSupport)
SolverParams & solverParams()
Get the solver parameters.
virtual std::size_t numNonlinearSystems() const override
Moose::LineSearchType _line_search
Definition: SolverParams.h:20
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
std::vector< std::pair< R1, R2 > > get(const std::string &param1, const std::string &param2) const
Combine two vector parameters into a single vector of pairs.
A struct for storing the various types of petsc options and values.
Definition: PetscSupport.h:38
Moose::MffdType _mffd_type
Definition: SolverParams.h:21
LineSearchType
Type of the line search.
Definition: MooseTypes.h:838
virtual unsigned int dimension() const
Returns MeshBase::mesh_dimension(), (not MeshBase::spatial_dimension()!) of the underlying libMesh me...
Definition: MooseMesh.C:2679
void processPetscPairs(const std::vector< std::pair< MooseEnumItem, std::string >> &petsc_pair_options, const unsigned int mesh_dimension, PetscOptions &petsc_options)
Populate name and value pairs in a given PetscOptions object using vectors of input arguments...
Definition: PetscSupport.C:643
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
Definition: MooseEnum.h:31
Moose::SolveType _type
Definition: SolverParams.h:19
NonlinearSystemBase & getNonlinearSystemBase(const unsigned int sys_num)
void processPetscFlags(const MultiMooseEnum &petsc_flags, PetscOptions &petsc_options)
Populate flags in a given PetscOptions object using a vector of input arguments.
Definition: PetscSupport.C:603
Class for containing MooseEnum item information.
Definition: MooseEnumItem.h:18
virtual System & system() override
Get the reference to the libMesh system.
IntRange< T > make_range(T beg, T end)
virtual MooseMesh & mesh() override
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type It sho...
means not set
Definition: MooseTypes.h:840
bool isParamValid(const std::string &name) const
This method returns parameters that have been initialized in one fashion or another, i.e.

◆ stringify() [1/2]

std::string Moose::PetscSupport::stringify ( const LineSearchType t)

Definition at line 95 of file PetscSupport.C.

Referenced by setSolverOptions().

96 {
97  switch (t)
98  {
99  case LS_BASIC:
100  return "basic";
101  case LS_DEFAULT:
102  return "default";
103  case LS_NONE:
104  return "none";
105  case LS_SHELL:
106  return "shell";
107  case LS_L2:
108  return "l2";
109  case LS_BT:
110  return "bt";
111  case LS_CP:
112  return "cp";
113  case LS_CONTACT:
114  return "contact";
115  case LS_PROJECT:
116  return "project";
117  case LS_INVALID:
118  mooseError("Invalid LineSearchType");
119  }
120  return "";
121 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
means not set
Definition: MooseTypes.h:840

◆ stringify() [2/2]

std::string Moose::PetscSupport::stringify ( const MffdType t)

Definition at line 124 of file PetscSupport.C.

125 {
126  switch (t)
127  {
128  case MFFD_WP:
129  return "wp";
130  case MFFD_DS:
131  return "ds";
132  case MFFD_INVALID:
133  mooseError("Invalid MffdType");
134  }
135  return "";
136 }
void mooseError(Args &&... args)
Emit an error message with the given stringified, concatenated args and terminate the application...
Definition: MooseError.h:299
means not set
Definition: MooseTypes.h:857