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 "ExplicitTimeIntegrator.h" 12 : #include "NonlinearSystem.h" 13 : #include "FEProblem.h" 14 : #include "PetscSupport.h" 15 : #include "KernelBase.h" 16 : #include "DGKernelBase.h" 17 : #include "ScalarKernelBase.h" 18 : #include "FVElementalKernel.h" 19 : #include "FVFluxKernel.h" 20 : #include "NodalKernelBase.h" 21 : 22 : // libMesh includes 23 : #include "libmesh/enum_convergence_flags.h" 24 : 25 : using namespace libMesh; 26 : 27 : InputParameters 28 43734 : ExplicitTimeIntegrator::validParams() 29 : { 30 43734 : InputParameters params = TimeIntegrator::validParams(); 31 : 32 43734 : MooseEnum solve_type("consistent lumped lump_preconditioned", "consistent"); 33 : 34 43734 : params.addParam<MooseEnum>( 35 : "solve_type", 36 : solve_type, 37 : "The way to solve the system. A 'consistent' solve uses the full mass matrix and actually " 38 : "needs to use a linear solver to solve the problem. 'lumped' uses a lumped mass matrix with " 39 : "a simple inversion - incredibly fast but may be less accurate. 'lump_preconditioned' uses " 40 : "the lumped mass matrix as a preconditioner for the 'consistent' solve"); 41 : 42 87468 : return params; 43 43734 : } 44 : 45 626 : ExplicitTimeIntegrator::ExplicitTimeIntegrator(const InputParameters & parameters) 46 : : TimeIntegrator(parameters), 47 : MeshChangedInterface(parameters), 48 : 49 626 : _solve_type(getParam<MooseEnum>("solve_type")), 50 626 : _explicit_residual(addVector("explicit_residual", false, PARALLEL)), 51 626 : _solution_update(addVector("solution_update", true, PARALLEL)), 52 1252 : _mass_matrix_diag_inverted(addVector("mass_matrix_diag_inverted", true, GHOSTED)) 53 : { 54 626 : _Ke_time_tag = _fe_problem.getMatrixTagID("TIME"); 55 : 56 : // This effectively changes the default solve_type to LINEAR instead of PJFNK, 57 : // so that it is valid to not supply solve_type in the Executioner block: 58 626 : if (_nl) 59 313 : _fe_problem.solverParams(_nl->number())._type = Moose::ST_LINEAR; 60 : 61 626 : if (_solve_type == LUMPED || _solve_type == LUMP_PRECONDITIONED) 62 156 : _ones = addVector("ones", true, PARALLEL); 63 : // don't set any of the common SNES-related petsc options to prevent unused option warnings 64 626 : Moose::PetscSupport::dontAddCommonSNESOptions(_fe_problem); 65 626 : } 66 : 67 : void 68 305 : ExplicitTimeIntegrator::initialSetup() 69 : { 70 305 : meshChanged(); 71 : 72 305 : if (_nl) 73 : { 74 305 : std::unordered_set<unsigned int> vars_to_check; 75 305 : if (!_vars.empty()) 76 0 : vars_to_check = _vars; 77 : else 78 622 : for (const auto i : make_range(_nl->nVariables())) 79 317 : vars_to_check.insert(i); 80 : 81 305 : const auto mass_matrix_tag_id = massMatrixTagID(); 82 305 : std::set<TagID> matrix_tags = {mass_matrix_tag_id}; 83 305 : auto fv_object_starting_query = _fe_problem.theWarehouse() 84 305 : .query() 85 305 : .template condition<AttribMatrixTags>(matrix_tags) 86 305 : .template condition<AttribSysNum>(_nl->number()); 87 : 88 618 : for (const auto var_id : vars_to_check) 89 : { 90 317 : const auto & var_name = _nl->system().variable_name(var_id); 91 317 : const bool field_var = _nl->hasVariable(var_name); 92 317 : if (!field_var) 93 : mooseAssert(_nl->hasScalarVariable(var_name), 94 : var_name << " should be either a field or scalar variable"); 95 317 : if (field_var) 96 : { 97 265 : const auto & field_var = _nl->getVariable(/*tid=*/0, var_name); 98 265 : if (field_var.isFV()) 99 : { 100 0 : std::vector<FVElementalKernel *> fv_elemental_kernels; 101 0 : auto var_query = fv_object_starting_query.clone().template condition<AttribVar>(var_id); 102 0 : auto var_query_clone = var_query.clone(); 103 0 : var_query.template condition<AttribSystem>("FVElementalKernel") 104 0 : .queryInto(fv_elemental_kernels); 105 0 : if (fv_elemental_kernels.size()) 106 0 : continue; 107 : 108 0 : std::vector<FVFluxKernel *> fv_flux_kernels; 109 0 : var_query_clone.template condition<AttribSystem>("FVFluxKernel") 110 0 : .queryInto(fv_flux_kernels); 111 0 : if (fv_flux_kernels.size()) 112 0 : continue; 113 0 : } 114 : else 115 : { 116 : // We are a finite element variable 117 265 : if (_nl->getKernelWarehouse() 118 265 : .getMatrixTagObjectWarehouse(mass_matrix_tag_id, 0) 119 265 : .hasVariableObjects(var_id)) 120 261 : continue; 121 4 : if (_nl->getDGKernelWarehouse() 122 4 : .getMatrixTagObjectWarehouse(mass_matrix_tag_id, 0) 123 4 : .hasVariableObjects(var_id)) 124 0 : continue; 125 4 : if (_nl->getNodalKernelWarehouse() 126 4 : .getMatrixTagObjectWarehouse(mass_matrix_tag_id, 0) 127 4 : .hasVariableObjects(var_id)) 128 0 : continue; 129 : } 130 : } 131 52 : else if (_nl->getScalarKernelWarehouse() 132 52 : .getMatrixTagObjectWarehouse(mass_matrix_tag_id, 0) 133 52 : .hasVariableObjects(var_id)) 134 52 : continue; 135 : 136 4 : mooseError("No objects contributing to the mass matrix were found for variable '", 137 : var_name, 138 : "'. Did you, e.g., forget a time derivative term?"); 139 : } 140 301 : } 141 301 : } 142 : 143 : void 144 279 : ExplicitTimeIntegrator::init() 145 : { 146 279 : if (_nl && _fe_problem.solverParams(_nl->number())._type != Moose::ST_LINEAR) 147 4 : mooseError( 148 : "The chosen time integrator requires 'solve_type = LINEAR' in the Executioner block."); 149 275 : } 150 : 151 : void 152 4334 : ExplicitTimeIntegrator::preSolve() 153 : { 154 4334 : } 155 : 156 : void 157 349 : ExplicitTimeIntegrator::meshChanged() 158 : { 159 : // Can only be done after the system is initialized 160 349 : if (_solve_type == LUMPED || _solve_type == LUMP_PRECONDITIONED) 161 78 : *_ones = 1.; 162 : 163 349 : if (_solve_type == CONSISTENT || _solve_type == LUMP_PRECONDITIONED) 164 300 : _linear_solver = LinearSolver<Number>::build(comm()); 165 : 166 349 : if (_solve_type == LUMP_PRECONDITIONED) 167 : { 168 29 : _preconditioner = std::make_unique<LumpedPreconditioner>(*_mass_matrix_diag_inverted); 169 29 : _linear_solver->attach_preconditioner(_preconditioner.get()); 170 29 : _linear_solver->init(); 171 : } 172 : 173 349 : if (_solve_type == CONSISTENT || _solve_type == LUMP_PRECONDITIONED) 174 300 : Moose::PetscSupport::setLinearSolverDefaults(_fe_problem, *_linear_solver); 175 349 : } 176 : 177 : bool 178 4466 : ExplicitTimeIntegrator::performExplicitSolve(SparseMatrix<Number> & mass_matrix) 179 : { 180 4466 : bool converged = false; 181 : 182 4466 : switch (_solve_type) 183 : { 184 1729 : case CONSISTENT: 185 : { 186 1729 : converged = solveLinearSystem(mass_matrix); 187 : 188 1729 : break; 189 : } 190 2443 : case LUMPED: 191 : { 192 : // Computes the sum of each row (lumping) 193 : // Note: This is actually how PETSc does it 194 : // It's not "perfectly optimal" - but it will be fast (and universal) 195 2443 : mass_matrix.vector_mult(*_mass_matrix_diag_inverted, *_ones); 196 : 197 : // "Invert" the diagonal mass matrix 198 2443 : _mass_matrix_diag_inverted->reciprocal(); 199 : 200 : // Multiply the inversion by the RHS 201 2443 : _solution_update->pointwise_mult(*_mass_matrix_diag_inverted, *_explicit_residual); 202 : 203 : // Check for convergence by seeing if there is a nan or inf 204 2443 : auto sum = _solution_update->sum(); 205 2443 : converged = std::isfinite(sum); 206 : 207 : // The linear iteration count remains zero 208 2443 : _n_linear_iterations = 0; 209 : 210 2443 : break; 211 : } 212 294 : case LUMP_PRECONDITIONED: 213 : { 214 294 : mass_matrix.vector_mult(*_mass_matrix_diag_inverted, *_ones); 215 294 : _mass_matrix_diag_inverted->reciprocal(); 216 : 217 294 : converged = solveLinearSystem(mass_matrix); 218 : 219 294 : break; 220 : } 221 0 : default: 222 0 : mooseError("Unknown solve_type in ExplicitTimeIntegrator."); 223 : } 224 : 225 4466 : return converged; 226 : } 227 : 228 : bool 229 2023 : ExplicitTimeIntegrator::solveLinearSystem(SparseMatrix<Number> & mass_matrix) 230 : { 231 2023 : auto & es = _fe_problem.es(); 232 : 233 : const auto num_its_and_final_tol = 234 2447 : _linear_solver->solve(mass_matrix, 235 2023 : *_solution_update, 236 2023 : *_explicit_residual, 237 212 : es.parameters.get<Real>("linear solver tolerance"), 238 212 : es.parameters.get<unsigned int>("linear solver maximum iterations")); 239 : 240 2023 : _n_linear_iterations += num_its_and_final_tol.first; 241 : 242 2023 : const bool converged = checkLinearConvergence(); 243 : 244 2023 : return converged; 245 : } 246 : 247 : bool 248 2023 : ExplicitTimeIntegrator::checkLinearConvergence() 249 : { 250 2023 : auto reason = _linear_solver->get_converged_reason(); 251 : 252 2023 : switch (reason) 253 : { 254 1879 : case CONVERGED_RTOL_NORMAL: 255 : case CONVERGED_ATOL_NORMAL: 256 : case CONVERGED_RTOL: 257 : case CONVERGED_ATOL: 258 : case CONVERGED_ITS: 259 : case CONVERGED_CG_NEG_CURVE: 260 : case CONVERGED_CG_CONSTRAINED: 261 : case CONVERGED_STEP_LENGTH: 262 : case CONVERGED_HAPPY_BREAKDOWN: 263 1879 : return true; 264 144 : case DIVERGED_NULL: 265 : case DIVERGED_ITS: 266 : case DIVERGED_DTOL: 267 : case DIVERGED_BREAKDOWN: 268 : case DIVERGED_BREAKDOWN_BICG: 269 : case DIVERGED_NONSYMMETRIC: 270 : case DIVERGED_INDEFINITE_PC: 271 : case DIVERGED_NAN: 272 : case DIVERGED_INDEFINITE_MAT: 273 : case CONVERGED_ITERATING: 274 : case DIVERGED_PCSETUP_FAILED: 275 144 : return false; 276 0 : default: 277 0 : mooseError("Unknown convergence reason in ExplicitTimeIntegrator."); 278 : } 279 : }