https://mooseframework.inl.gov
SlepcSupport.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "MultiMooseEnum.h"
13 #include "libmesh/libmesh_config.h"
14 
15 #ifdef LIBMESH_HAVE_SLEPC
16 
17 #include "libmesh/petsc_macro.h"
18 
19 #include "Moose.h"
20 #include "SystemBase.h"
21 
22 /* We need this in order to implement moose PC */
23 #include <petsc/private/pcimpl.h>
24 #include <slepceps.h>
25 #include <slepc/private/epsimpl.h>
26 /* In order to use libMesh preconditioner */
27 #include "libmesh/linear_solver.h"
28 #include "libmesh/preconditioner.h"
29 
30 class EigenProblem;
31 class InputParameters;
32 class SolverParams;
33 class MultiMooseEnum;
34 
35 namespace Moose
36 {
37 namespace SlepcSupport
38 {
46 
51 
55 void storeSolveType(FEProblemBase & fe_problem, const InputParameters & params);
56 
60 void setEigenProblemSolverParams(EigenProblem & eigen_problem, const InputParameters & params);
61 
66 void slepcSetOptions(EigenProblem & eigen_problem,
67  SolverParams & solver_params,
68  const InputParameters & params);
69 
73 void setSlepcEigenSolverTolerances(EigenProblem & eigen_problem,
74  const SolverParams & solver_params,
75  const InputParameters & params);
76 
80 void setFreeNonlinearPowerIterations(unsigned int free_power_iterations);
81 
82 /*
83  * Set SLEPc/PETSc options to turn the eigen-solver back to a regular Newton solver
84  */
86 
103  SNES snes, Vec x, Mat eigen_mat, SparseMatrix<Number> & all_dofs_mat, void * ctx, TagID tag);
104 
109 void moosePetscSNESFormMatricesTags(SNES snes,
110  Vec x,
111  std::vector<Mat> & eigen_mats,
112  std::vector<SparseMatrix<Number> *> & all_dofs_mats,
113  void * ctx,
114  const std::set<TagID> & tags);
115 
119 PetscErrorCode mooseSlepcEigenFormJacobianA(SNES snes, Vec x, Mat jac, Mat pc, void * ctx);
120 
124 PetscErrorCode mooseSlepcEigenFormJacobianB(SNES snes, Vec x, Mat jac, Mat pc, void * ctx);
125 
129 PetscErrorCode mooseSlepcEigenFormFunctionA(SNES snes, Vec x, Vec r, void * ctx);
130 
134 PetscErrorCode mooseSlepcEigenFormFunctionB(SNES snes, Vec x, Vec r, void * ctx);
135 
139 PetscErrorCode mooseSlepcEigenFormFunctionAB(SNES snes, Vec x, Vec Ax, Vec Bx, void * ctx);
140 
145 PetscErrorCode mooseSlepcStoppingTest(EPS eps,
146  PetscInt its,
147  PetscInt max_it,
148  PetscInt nconv,
149  PetscInt nev,
150  EPSConvergedReason * reason,
151  void * ctx);
152 
156 PetscErrorCode mooseSlepcEPSGetSNES(EPS eps, SNES * snes);
157 
161 PetscErrorCode mooseSlepcEPSMonitor(EPS eps,
162  PetscInt its,
163  PetscInt nconv,
164  PetscScalar * eigr,
165  PetscScalar * eigi,
166  PetscReal * errest,
167  PetscInt nest,
168  void * mctx);
169 
173 PetscErrorCode mooseSlepcEPSSNESSetUpOptionPrefix(EPS eps);
174 
179 PetscErrorCode mooseSlepcEPSSNESSetCustomizePC(EPS eps);
180 
185 PetscErrorCode mooseSlepcEPSSNESKSPSetPCSide(FEProblemBase & problem, EPS eps);
186 
190 void attachCallbacksToMat(EigenProblem & eigen_problem, Mat mat, bool eigen);
191 
195 PetscErrorCode mooseMatMult_Eigen(Mat mat, Vec x, Vec y);
196 
200 PetscErrorCode mooseMatMult_NonEigen(Mat mat, Vec x, Vec y);
201 
205 void setOperationsForShellMat(EigenProblem & eigen_problem, Mat mat, bool eigen);
206 
210 PETSC_EXTERN PetscErrorCode PCCreate_MoosePC(PC pc);
211 
215 PETSC_EXTERN PetscErrorCode registerPCToPETSc();
216 
220 PetscErrorCode PCDestroy_MoosePC(PC pc);
221 
225 PetscErrorCode PCView_MoosePC(PC pc, PetscViewer viewer);
226 
230 PetscErrorCode PCApply_MoosePC(PC pc, Vec x, Vec y);
231 
235 PetscErrorCode PCSetUp_MoosePC(PC pc);
236 
240 PetscErrorCode mooseSlepcEigenFormFunctionMFFD(void * ctx, Vec x, Vec r);
241 
242 } // namespace SlepcSupport
243 } // namespace moose
244 
245 #endif // LIBMESH_HAVE_SLEPC
int eps(unsigned int i, unsigned int j)
2D version
PetscErrorCode mooseSlepcEigenFormFunctionMFFD(void *ctx, Vec x, Vec r)
Function call for MFFD.
Definition: SlepcSupport.C:749
unsigned int TagID
Definition: MooseTypes.h:210
PetscErrorCode mooseSlepcEPSSNESSetCustomizePC(EPS eps)
Attach a customized PC.
PetscErrorCode mooseSlepcStoppingTest(EPS eps, PetscInt its, PetscInt max_it, PetscInt nconv, PetscInt nev, EPSConvergedReason *reason, void *ctx)
A customized convergence checker.
PetscErrorCode mooseMatMult_Eigen(Mat mat, Vec x, Vec y)
Implement MatMult via function evaluation for Bx.
PetscErrorCode mooseSlepcEigenFormJacobianB(SNES snes, Vec x, Mat jac, Mat pc, void *ctx)
Form Jacobian matrix B.
Definition: SlepcSupport.C:861
void clearFreeNonlinearPowerIterations(const InputParameters &params)
Definition: SlepcSupport.C:419
PetscErrorCode mooseSlepcEigenFormFunctionAB(SNES snes, Vec x, Vec Ax, Vec Bx, void *ctx)
Form function residual Ax-Bx.
Definition: SlepcSupport.C:985
The main MOOSE class responsible for handling user-defined parameters in almost every MOOSE system...
PetscErrorCode PCApply_MoosePC(PC pc, Vec x, Vec y)
Preconditioner application.
PetscErrorCode PCDestroy_MoosePC(PC pc)
Destroy preconditioner.
PetscErrorCode mooseSlepcEigenFormFunctionA(SNES snes, Vec x, Vec r, void *ctx)
Form function residual Ax.
Definition: SlepcSupport.C:911
Specialization of SubProblem for solving nonlinear equations plus auxiliary equations.
PetscErrorCode mooseMatMult_NonEigen(Mat mat, Vec x, Vec y)
Implement MatMult via function evaluation for Ax.
void setSlepcEigenSolverTolerances(EigenProblem &eigen_problem, const SolverParams &solver_params, const InputParameters &params)
Control eigen solver tolerances via SLEPc options.
Definition: SlepcSupport.C:129
InputParameters getSlepcValidParams(InputParameters &params)
Definition: SlepcSupport.C:38
PetscErrorCode PCSetUp_MoosePC(PC pc)
Setup preconditioner.
void attachCallbacksToMat(EigenProblem &eigen_problem, Mat mat, bool eigen)
Attach call backs to mat.
void slepcSetOptions(EigenProblem &eigen_problem, SolverParams &solver_params, const InputParameters &params)
Push all SLEPc/PETSc options into SLEPc/PETSc side.
Definition: SlepcSupport.C:550
PETSC_EXTERN PetscErrorCode PCCreate_MoosePC(PC pc)
Create a preconditioner from moose side.
PetscErrorCode PCView_MoosePC(PC pc, PetscViewer viewer)
View preconditioner.
PetscErrorCode mooseSlepcEigenFormJacobianA(SNES snes, Vec x, Mat jac, Mat pc, void *ctx)
Form Jacobian matrix A.
Definition: SlepcSupport.C:775
void moosePetscSNESFormMatricesTags(SNES snes, Vec x, std::vector< Mat > &eigen_mats, std::vector< SparseMatrix< Number > *> &all_dofs_mats, void *ctx, const std::set< TagID > &tags)
Form multiple matrices for multiple tags.
PetscErrorCode mooseSlepcEigenFormFunctionB(SNES snes, Vec x, Vec r, void *ctx)
Form function residual Bx.
Definition: SlepcSupport.C:943
PetscErrorCode mooseSlepcEPSMonitor(EPS eps, PetscInt its, PetscInt nconv, PetscScalar *eigr, PetscScalar *eigi, PetscReal *errest, PetscInt nest, void *mctx)
A customized solver monitor to print out eigenvalue.
PetscErrorCode mooseSlepcEPSGetSNES(EPS eps, SNES *snes)
Retrieve SNES from EPS.
void setFreeNonlinearPowerIterations(unsigned int free_power_iterations)
Set SLEPc/PETSc options to trigger free power iteration.
Definition: SlepcSupport.C:404
PetscErrorCode mooseSlepcEPSSNESKSPSetPCSide(FEProblemBase &problem, EPS eps)
Allow users to specify PC side.
PetscErrorCode mooseSlepcEPSSNESSetUpOptionPrefix(EPS eps)
Get rid of prefix "-eps_power" for SNES, KSP, PC, etc.
void moosePetscSNESFormMatrixTag(SNES snes, Vec x, Mat eigen_mat, SparseMatrix< Number > &all_dofs_mat, void *ctx, TagID tag)
Form matrix according to tag.
void * ctx
void setOperationsForShellMat(EigenProblem &eigen_problem, Mat mat, bool eigen)
Set operations to shell mat.
void setEigenProblemSolverParams(EigenProblem &eigen_problem, const InputParameters &params)
Retrieve eigen problem params from &#39;params&#39;, and then set these params into SolverParams.
Definition: SlepcSupport.C:212
MOOSE now contains C++17 code, so give a reasonable error message stating what the user can do to add...
This is a "smart" enum class intended to replace many of the shortcomings in the C++ enum type...
Problem for solving eigenvalue problems.
Definition: EigenProblem.h:21
void storeSolveType(FEProblemBase &fe_problem, const InputParameters &params)
Set solve type into eigen problem (solverParams)
Definition: SlepcSupport.C:286
PETSC_EXTERN PetscErrorCode registerPCToPETSc()
Let PETSc know there is a preconditioner.
InputParameters getSlepcEigenProblemValidParams()
Retrieve valid params that allow users to specify eigen problem configuration.
Definition: SlepcSupport.C:67