Line data Source code
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 : // MOOSE includes 11 : #include "EigenProblem.h" 12 : #include "Factory.h" 13 : #include "MooseApp.h" 14 : #include "NonlinearEigenSystem.h" 15 : #include "PetscSupport.h" 16 : #include "SlepcSupport.h" 17 : #include "UserObject.h" 18 : 19 : #include "libmesh/petsc_solver_exception.h" 20 : 21 : // Needed for LIBMESH_CHECK_ERR 22 : using libMesh::PetscSolverException; 23 : 24 : InputParameters 25 4188 : EigenProblemSolve::validParams() 26 : { 27 4188 : InputParameters params = Executioner::validParams(); 28 : 29 4188 : params.addClassDescription( 30 : "Solve object for a standard/generalized linear or nonlinear eigenvalue problem"); 31 : 32 4188 : params += FEProblemSolve::validParams(); 33 : 34 12564 : params.addParam<bool>( 35 : "matrix_free", 36 8376 : false, 37 : "Whether or not to use a matrix free fashion to form operators. " 38 : "If true, shell matrices will be used and meanwhile a preconditioning matrix" 39 : "may be formed as well."); 40 : 41 12564 : params.addParam<bool>( 42 : "precond_matrix_free", 43 8376 : false, 44 : "Whether or not to use a matrix free fashion for forming the preconditioning matrix. " 45 : "If true, a shell matrix will be used for preconditioner."); 46 : 47 12564 : params.addParam<bool>("constant_matrices", 48 8376 : false, 49 : "Whether or not to use constant matrices so that we can use them to form " 50 : "residuals on both linear and " 51 : "nonlinear iterations"); 52 : 53 8376 : params.addParam<bool>( 54 : "precond_matrix_includes_eigen", 55 8376 : false, 56 : "Whether or not to include eigen kernels in the preconditioning matrix. " 57 : "If true, the preconditioning matrix could be singular with the converged eigenvalue if the " 58 : "full matrix is assembled and the derivative of eigenvalue with respect to the solution " 59 : "vector is not considered."); 60 : 61 12564 : params.addPrivateParam<bool>("_use_eigen_value", true); 62 : 63 16752 : params.addParam<Real>("initial_eigenvalue", 1, "Initial eigenvalue"); 64 16752 : params.addParam<PostprocessorName>( 65 : "normalization", "Postprocessor evaluating norm of eigenvector for normalization"); 66 16752 : params.addParam<Real>("normal_factor", 67 : "Normalize eigenvector to make a defined norm equal to this factor"); 68 : 69 12564 : params.addParam<bool>("auto_initialization", 70 8376 : true, 71 : "If true, we will set an initial eigen vector in moose, otherwise EPS " 72 : "solver will initialize eigen vector"); 73 : 74 16752 : params.addParamNamesToGroup("matrix_free precond_matrix_free constant_matrices " 75 : "precond_matrix_includes_eigen", 76 : "Matrix and Matrix-Free"); 77 16752 : params.addParamNamesToGroup("initial_eigenvalue auto_initialization", 78 : "Eigenvector and eigenvalue initialization"); 79 12564 : params.addParamNamesToGroup("normalization normal_factor", "Solution normalization"); 80 : 81 : // If Newton and Inverse Power is combined in SLEPc side 82 8376 : params.addPrivateParam<bool>("_newton_inverse_power", false); 83 : 84 : // Add slepc options and eigen problems 85 : #ifdef LIBMESH_HAVE_SLEPC 86 4188 : Moose::SlepcSupport::getSlepcValidParams(params); 87 : 88 4188 : params += Moose::SlepcSupport::getSlepcEigenProblemValidParams(); 89 : #endif 90 4188 : return params; 91 0 : } 92 : 93 565 : EigenProblemSolve::EigenProblemSolve(Executioner & ex) 94 : : FEProblemSolve(ex), 95 1695 : _eigen_problem(*getCheckedPointerParam<EigenProblem *>( 96 : "_eigen_problem", "This might happen if you don't have a mesh")), 97 1148 : _normalization(isParamValid("normalization") ? &getPostprocessorValue("normalization") 98 565 : : nullptr) 99 : { 100 : // Extract and store SLEPc options 101 : #ifdef LIBMESH_HAVE_SLEPC 102 : mooseAssert(_problem.numSolverSystems() == 1, 103 : "The Eigenvalue executioner only currently supports a single solver system."); 104 : 105 565 : const auto & params = ex.parameters(); 106 565 : Moose::SlepcSupport::storeSolveType(_eigen_problem, params); 107 : 108 565 : Moose::SlepcSupport::setEigenProblemSolverParams(_eigen_problem, params); 109 1124 : _eigen_problem.setEigenproblemType( 110 562 : _eigen_problem.solverParams(/*solver_sys_num=*/0)._eigen_problem_type); 111 : 112 : // pass two control parameters to eigen problem 113 1124 : _eigen_problem.solverParams(/*solver_sys_num=*/0)._free_power_iterations = 114 1124 : getParam<unsigned int>("free_power_iterations"); 115 1124 : _eigen_problem.solverParams(/*solver_sys_num=*/0)._extra_power_iterations = 116 1124 : getParam<unsigned int>("extra_power_iterations"); 117 : 118 2792 : if (!isParamValid("normalization") && isParamValid("normal_factor")) 119 0 : paramError("normal_factor", 120 : "Cannot set scaling factor without defining normalization postprocessor."); 121 : 122 1686 : if (isParamValid("normalization")) 123 : { 124 18 : const auto & normpp = getParam<PostprocessorName>("normalization"); 125 27 : if (isParamValid("normal_factor")) 126 27 : _eigen_problem.setNormalization(normpp, getParam<Real>("normal_factor")); 127 : else 128 0 : _eigen_problem.setNormalization(normpp); 129 : } 130 : 131 1124 : _eigen_problem.setInitialEigenvalue(getParam<Real>("initial_eigenvalue")); 132 : 133 : // Set a flag to nonlinear eigen system 134 562 : _eigen_problem.getNonlinearEigenSystem(/*nl_sys_num=*/0) 135 1124 : .precondMatrixIncludesEigenKernels(getParam<bool>("precond_matrix_includes_eigen")); 136 : #else 137 : mooseError("SLEPc is required to use Eigenvalue executioner, please use '--download-slepc in " 138 : "PETSc configuration'"); 139 : #endif 140 : // SLEPc older than 3.13.0 can not take initial guess from moose 141 : // It may generate converge issues 142 : #if PETSC_RELEASE_LESS_THAN(3, 13, 0) 143 : mooseError( 144 : "Please use SLEPc-3.13.0 or higher. Old versions of SLEPc likely produce bad convergence"); 145 : #endif 146 : 147 : // To avoid petsc unused option warnings, ensure we do not set irrelevant options. 148 562 : Moose::PetscSupport::dontAddLinearConvergedReason(_problem); 149 562 : Moose::PetscSupport::dontAddNonlinearConvergedReason(_problem); 150 1124 : Moose::PetscSupport::dontAddPetscFlag("-mat_mffd_type", _problem.getPetscOptions()); 151 1124 : Moose::PetscSupport::dontAddPetscFlag("-st_ksp_atol", _problem.getPetscOptions()); 152 1124 : Moose::PetscSupport::dontAddPetscFlag("-st_ksp_max_it", _problem.getPetscOptions()); 153 1124 : Moose::PetscSupport::dontAddPetscFlag("-st_ksp_rtol", _problem.getPetscOptions()); 154 562 : if (!_eigen_problem.solverParams()._eigen_matrix_free) 155 243 : Moose::PetscSupport::dontAddPetscFlag("-snes_mf_operator", _problem.getPetscOptions()); 156 562 : } 157 : 158 : void 159 550 : EigenProblemSolve::initialSetup() 160 : { 161 1650 : if (isParamValid("normalization")) 162 : { 163 18 : const auto & normpp = getParam<PostprocessorName>("normalization"); 164 9 : const auto & exec = _eigen_problem.getUserObject<UserObject>(normpp).getExecuteOnEnum(); 165 9 : if (!exec.isValueSet(EXEC_LINEAR)) 166 0 : mooseError("Normalization postprocessor ", normpp, " requires execute_on = 'linear'"); 167 : } 168 : 169 : #ifdef LIBMESH_HAVE_SLEPC 170 : // Options need to be setup once only 171 550 : if (!_eigen_problem.petscOptionsInserted()) 172 : { 173 : // Parent application has the default data base 174 541 : if (!_app.isUltimateMaster()) 175 24 : LibmeshPetscCall(PetscOptionsPush(_eigen_problem.petscOptionsDatabase())); 176 541 : Moose::SlepcSupport::slepcSetOptions( 177 541 : _eigen_problem, _eigen_problem.solverParams(/*eigen_sys_num=*/0), _pars); 178 541 : if (!_app.isUltimateMaster()) 179 24 : LibmeshPetscCall(PetscOptionsPop()); 180 541 : _eigen_problem.petscOptionsInserted() = true; 181 : } 182 : #endif 183 550 : }