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 : #pragma once 11 : 12 : #include "TimeIntegrator.h" 13 : #include "MeshChangedInterface.h" 14 : 15 : #include "libmesh/linear_solver.h" 16 : #include "libmesh/preconditioner.h" 17 : #include "libmesh/numeric_vector.h" 18 : 19 : #include "LumpedPreconditioner.h" 20 : 21 : /** 22 : * Base class for explicit time integrators that are implemented without using 23 : * a nonlinear solver. 24 : */ 25 : class ExplicitTimeIntegrator : public TimeIntegrator, public MeshChangedInterface 26 : { 27 : public: 28 : static InputParameters validParams(); 29 : 30 : ExplicitTimeIntegrator(const InputParameters & parameters); 31 : 32 : virtual void initialSetup() override; 33 : virtual void init() override; 34 : virtual void preSolve() override; 35 : virtual void meshChanged() override; 36 0 : virtual bool isExplicit() const override { return true; } 37 : 38 : protected: 39 : enum SolveType 40 : { 41 : CONSISTENT, 42 : LUMPED, 43 : LUMP_PRECONDITIONED 44 : }; 45 : 46 : /** 47 : * Solves a linear system using the chosen solve type 48 : * 49 : * @param[in] mass_matrix Mass matrix 50 : */ 51 : virtual bool performExplicitSolve(SparseMatrix<Number> & mass_matrix); 52 : 53 : /** 54 : * Solves a linear system 55 : * 56 : * @param[in] mass_matrix Mass matrix 57 : */ 58 : bool solveLinearSystem(SparseMatrix<Number> & mass_matrix); 59 : 60 : /** 61 : * Check for the linear solver convergence 62 : */ 63 : bool checkLinearConvergence(); 64 : 65 : /** 66 : * @returns the mass matrix tag ID 67 : */ 68 305 : virtual TagID massMatrixTagID() const { return _Ke_time_tag; } 69 : 70 : /// Solve type for how mass matrix is handled 71 : MooseEnum _solve_type; 72 : 73 : /// Residual used for the RHS 74 : NumericVector<Real> * _explicit_residual; 75 : 76 : /// Solution vector for the linear solve 77 : NumericVector<Real> * _solution_update; 78 : 79 : /// Diagonal of the lumped mass matrix (and its inversion) 80 : NumericVector<Real> * _mass_matrix_diag_inverted; 81 : 82 : /// Vector of 1's to help with creating the lumped mass matrix 83 : NumericVector<Real> * _ones; 84 : 85 : /// For computing the mass matrix 86 : TagID _Ke_time_tag; 87 : 88 : /// For solving with the consistent matrix 89 : std::unique_ptr<libMesh::LinearSolver<Number>> _linear_solver; 90 : 91 : /// For solving with lumped preconditioning 92 : std::unique_ptr<LumpedPreconditioner> _preconditioner; 93 : 94 : /// Save off current time to reset it back and forth 95 : Real _current_time; 96 : };