Line data Source code
1 : // The libMesh Finite Element Library. 2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 3 : 4 : // This library is free software; you can redistribute it and/or 5 : // modify it under the terms of the GNU Lesser General Public 6 : // License as published by the Free Software Foundation; either 7 : // version 2.1 of the License, or (at your option) any later version. 8 : 9 : // This library is distributed in the hope that it will be useful, 10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of 11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 12 : // Lesser General Public License for more details. 13 : 14 : // You should have received a copy of the GNU Lesser General Public 15 : // License along with this library; if not, write to the Free Software 16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 17 : 18 : 19 : 20 : #ifndef LIBMESH_LINEAR_IMPLICIT_SYSTEM_H 21 : #define LIBMESH_LINEAR_IMPLICIT_SYSTEM_H 22 : 23 : // Local Includes 24 : #include "libmesh/implicit_system.h" 25 : 26 : // C++ includes 27 : #include <cstddef> 28 : 29 : namespace libMesh 30 : { 31 : 32 : 33 : // Forward Declarations 34 : template <typename T> class LinearSolver; 35 : template <typename T> class ShellMatrix; 36 : 37 : 38 : /** 39 : * \brief Manages consistently variables, degrees of freedom, coefficient 40 : * vectors, matrices and linear solvers for implicit systems. 41 : * 42 : * An implicit system is a system that requires the solution of a 43 : * system of equations. This class has the ability to create and use a 44 : * linear solver to solve the system. 45 : * 46 : * The matrix LinearImplicitSystem::matrix and the vector 47 : * LinearImplicitSystem::rhs should be filled during assembly. 48 : * 49 : * \note Additional vectors/matrices can be added via parent class 50 : * interfaces. 51 : * 52 : * \author Benjamin Kirk 53 : * \date 2005 54 : */ 55 1240 : class LinearImplicitSystem : public ImplicitSystem 56 : { 57 : public: 58 : 59 : /** 60 : * Constructor. 61 : */ 62 : LinearImplicitSystem (EquationSystems & es, 63 : const std::string & name, 64 : const unsigned int number); 65 : 66 : /** 67 : * Special functions. 68 : * - This class has the same restrictions/defaults as its base class. 69 : * - The destructor is defaulted out-of-line. 70 : */ 71 : LinearImplicitSystem (const LinearImplicitSystem &) = delete; 72 : LinearImplicitSystem & operator= (const LinearImplicitSystem &) = delete; 73 : LinearImplicitSystem (LinearImplicitSystem &&) = default; 74 : LinearImplicitSystem & operator= (LinearImplicitSystem &&) = delete; 75 : virtual ~LinearImplicitSystem (); 76 : 77 : /** 78 : * The type of system. 79 : */ 80 : typedef LinearImplicitSystem sys_type; 81 : 82 : /** 83 : * The type of the parent. 84 : */ 85 : typedef ImplicitSystem Parent; 86 : 87 : /** 88 : * \returns A reference to *this. 89 : */ 90 : sys_type & system () { return *this; } 91 : 92 : /** 93 : * Clear all the data structures associated with 94 : * the system. 95 : */ 96 : virtual void clear () override; 97 : 98 : /** 99 : * Initializes new data members of the system 100 : */ 101 : virtual void init_data () override; 102 : 103 : /** 104 : * Reinitializes the member data fields associated with 105 : * the system, so that, e.g., \p assemble() may be used. 106 : */ 107 : virtual void reinit () override; 108 : 109 : /** 110 : * Prepares \p matrix and \p _dof_map for matrix assembly. 111 : * Does not actually assemble anything. For matrix assembly, 112 : * use the \p assemble() in derived classes. 113 : * Should be overridden in derived classes. 114 : */ 115 38939 : virtual void assemble () override { ImplicitSystem::assemble(); } 116 : 117 : /** 118 : * After calling this method, any solve will be limited to the given 119 : * subset. To disable this mode, call this method with \p subset 120 : * being a \p nullptr. 121 : */ 122 : virtual void restrict_solve_to (const SystemSubset * subset, 123 : const SubsetSolveMode subset_solve_mode=SUBSET_ZERO) override; 124 : 125 : /** 126 : * Assembles & solves the linear system A*x=b. 127 : */ 128 : virtual void solve () override; 129 : 130 : /** 131 : * \returns A pointer to a linear solver appropriate for use in 132 : * adjoint and/or sensitivity solves 133 : */ 134 : virtual LinearSolver<Number> * get_linear_solver() const override; 135 : 136 : /** 137 : * Assembles a residual in \p rhs and/or a jacobian in \p matrix, 138 : * as requested. 139 : */ 140 : virtual void assembly(bool get_residual, 141 : bool get_jacobian, 142 : bool apply_heterogeneous_constraints = false, 143 : bool apply_no_constraints = false) override; 144 : 145 : /** 146 : * \returns \p "LinearImplicit". Helps in identifying 147 : * the system type in an equation system file. 148 : */ 149 10754 : virtual std::string system_type () const override { return "LinearImplicit"; } 150 : 151 : /** 152 : * \returns The number of iterations 153 : * taken for the most recent linear solve. 154 : */ 155 0 : unsigned int n_linear_iterations() const { return _n_linear_iterations; } 156 : 157 : /** 158 : * \returns The final residual for the linear system solve. 159 : */ 160 0 : Real final_linear_residual() const { return _final_linear_residual; } 161 : 162 : /** 163 : * This function enables the user to provide a shell matrix, i.e. a 164 : * matrix that is not stored element-wise, but as a function. When 165 : * you register your shell matrix using this function, calling \p 166 : * solve() will no longer use the \p matrix member but the 167 : * registered shell matrix instead. You can reset this behaviour to 168 : * its original state by supplying a \p nullptr to this 169 : * function. 170 : */ 171 : void attach_shell_matrix (ShellMatrix<Number> * shell_matrix); 172 : 173 : /** 174 : * Detaches a shell matrix. Same as \p attach_shell_matrix(nullptr). 175 : */ 176 : void detach_shell_matrix () { attach_shell_matrix(nullptr); } 177 : 178 : /** 179 : * \returns A pointer to the currently attached shell matrix, if any, 180 : * otherwise \p nullptr. 181 : */ 182 : ShellMatrix<Number> * get_shell_matrix() { return _shell_matrix; } 183 : 184 : virtual void create_static_condensation() override; 185 : 186 : protected: 187 : 188 : /** 189 : * The number of linear iterations required to solve the linear 190 : * system Ax=b. 191 : */ 192 : unsigned int _n_linear_iterations; 193 : 194 : /** 195 : * The final residual for the linear system Ax=b. 196 : */ 197 : Real _final_linear_residual; 198 : 199 : /** 200 : * User supplies shell matrix or \p nullptr if no shell matrix is used. 201 : */ 202 : ShellMatrix<Number> * _shell_matrix; 203 : 204 : /** 205 : * The current subset on which to solve (or \p nullptr if none). 206 : */ 207 : const SystemSubset * _subset; 208 : 209 : /** 210 : * If restrict-solve-to-subset mode is active, this member decides 211 : * what happens with the dofs outside the subset. 212 : */ 213 : SubsetSolveMode _subset_solve_mode; 214 : }; 215 : 216 : } // namespace libMesh 217 : 218 : #endif // LIBMESH_LINEAR_IMPLICIT_SYSTEM_H