Loading [MathJax]/extensions/tex2jax.js
libMesh
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
eigen_time_solver.C
Go to the documentation of this file.
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 #include "libmesh/libmesh_config.h"
21 #ifdef LIBMESH_HAVE_SLEPC
22 
23 #include "libmesh/diff_system.h"
24 #include "libmesh/eigen_time_solver.h"
25 #include "libmesh/eigen_solver.h"
26 #include "libmesh/sparse_matrix.h"
27 #include "libmesh/enum_eigen_solver_type.h"
28 
29 namespace libMesh
30 {
31 
32 // Constructor
34  : Parent(s),
35  eigen_solver (EigenSolver<Number>::build(s.comm())),
36  tol(pow(TOLERANCE, 5./3.)),
37  maxits(1000),
38  n_eigenpairs_to_compute(5),
39  n_basis_vectors_to_use(3*n_eigenpairs_to_compute),
40  n_converged_eigenpairs(0),
41  n_iterations_reqd(0)
42 {
43  libmesh_experimental();
44  eigen_solver->set_eigenproblem_type(GHEP);//or GNHEP
45  eigen_solver->set_position_of_spectrum(LARGEST_MAGNITUDE);
46 }
47 
49 
51 {
52  // empty...
53 }
54 
56 {
57  // Add matrix "B" to _system if not already there.
58  // The user may have already added a matrix "B" before
59  // calling the System initialization. This would be
60  // necessary if e.g. the System originally started life
61  // with a different type of TimeSolver and only later
62  // had its TimeSolver changed to an EigenTimeSolver.
63  if (!_system.have_matrix("B"))
64  _system.add_matrix("B");
65 }
66 
68 {
69  // The standard implementation is basically to call:
70  // _diff_solver->solve();
71  // which internally assembles (when necessary) matrices and vectors
72  // and calls linear solver software while also doing Newton steps (see newton_solver.C)
73  //
74  // The element_residual and side_residual functions below control
75  // what happens in the interior of the element assembly loops.
76  // We have a system reference, so it's possible to call _system.assembly()
77  // ourselves if we want to...
78  //
79  // Interestingly, for the EigenSolver we don't need residuals...just Jacobians.
80  // The Jacobian should therefore always be requested, and always return
81  // jacobian_computed as being true.
82 
83  // The basic plan of attack is:
84  // .) Construct the Jacobian using _system.assembly(true,true) as we
85  // would for a steady system. Use a flag in this class to
86  // control behavior in element_residual and side_residual
87  // .) Swap _system.matrix to matrix "B" (be sure to add this extra matrix during init)
88  // .) Call _system.assembly(true,true) again, using the flag in element_residual
89  // and side_residual to only get the mass matrix terms.
90  // .) Send A and B to Steffen's EigenSolver interface.
91 
92  // Assemble the spatial part (matrix A) of the operator
93  if (!this->quiet)
94  libMesh::out << "Assembling matrix A." << std::endl;
96  this->now_assembling = Matrix_A;
97  _system.assembly(true, true);
98  //_system.matrix->print_matlab("matrix_A.m");
99 
100  // Point the system's matrix at B, call assembly again.
101  if (!this->quiet)
102  libMesh::out << "Assembling matrix B." << std::endl;
103  _system.matrix = &( _system.get_matrix ("B") );
104  this->now_assembling = Matrix_B;
105  _system.assembly(true, true);
106  //_system.matrix->print_matlab("matrix_B.m");
107 
108  // Send matrices A, B to Steffen's SlepcEigenSolver interface
109  //libmesh_here();
110  if (!this->quiet)
111  libMesh::out << "Calling the EigenSolver." << std::endl;
112  std::tie(this->n_converged_eigenpairs, this->n_iterations_reqd) =
113  eigen_solver->solve_generalized (_system.get_system_matrix(),
114  _system.get_matrix ("B"),
117  tol,
118  maxits);
119 }
120 
121 
122 
123 bool EigenTimeSolver::element_residual(bool request_jacobian,
124  DiffContext & context)
125 {
126  // The EigenTimeSolver always computes jacobians!
127  libmesh_assert (request_jacobian);
128 
129  context.elem_solution_rate_derivative = 1;
130  context.elem_solution_derivative = 1;
131 
132  // Assemble the operator for the spatial part.
133  if (now_assembling == Matrix_A)
134  {
135  bool jacobian_computed =
136  _system.get_physics()->element_time_derivative(request_jacobian, context);
137 
138  // The user shouldn't compute a jacobian unless requested
139  libmesh_assert(request_jacobian || !jacobian_computed);
140 
141  bool jacobian_computed2 =
142  _system.get_physics()->element_constraint(jacobian_computed, context);
143 
144  // The user shouldn't compute a jacobian unless requested
145  libmesh_assert (jacobian_computed || !jacobian_computed2);
146 
147  return jacobian_computed && jacobian_computed2;
148 
149  }
150 
151  // Assemble the mass matrix operator
152  else if (now_assembling == Matrix_B)
153  {
154  bool mass_jacobian_computed =
155  _system.get_physics()->mass_residual(request_jacobian, context);
156 
157  // Scale Jacobian by -1 to get positive matrix from new negative
158  // mass_residual convention
159  context.get_elem_jacobian() *= -1.0;
160 
161  return mass_jacobian_computed;
162  }
163 
164  else
165  libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
166 }
167 
168 
169 
170 bool EigenTimeSolver::side_residual(bool request_jacobian,
171  DiffContext & context)
172 {
173  // The EigenTimeSolver always requests jacobians?
174  //libmesh_assert (request_jacobian);
175 
176  context.elem_solution_rate_derivative = 1;
177  context.elem_solution_derivative = 1;
178 
179  // Assemble the operator for the spatial part.
180  if (now_assembling == Matrix_A)
181  {
182  bool jacobian_computed =
183  _system.get_physics()->side_time_derivative(request_jacobian, context);
184 
185  // The user shouldn't compute a jacobian unless requested
186  libmesh_assert (request_jacobian || !jacobian_computed);
187 
188  bool jacobian_computed2 =
189  _system.get_physics()->side_constraint(jacobian_computed, context);
190 
191  // The user shouldn't compute a jacobian unless requested
192  libmesh_assert (jacobian_computed || !jacobian_computed2);
193 
194  return jacobian_computed && jacobian_computed2;
195 
196  }
197 
198  // There is now a "side" equivalent for the mass matrix
199  else if (now_assembling == Matrix_B)
200  {
201  bool mass_jacobian_computed =
202  _system.get_physics()->side_mass_residual(request_jacobian, context);
203 
204  // Scale Jacobian by -1 to get positive matrix from new negative
205  // mass_residual convention
206  context.get_elem_jacobian() *= -1.0;
207 
208  return mass_jacobian_computed;
209  }
210 
211  else
212  libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
213 }
214 
215 
216 
217 bool EigenTimeSolver::nonlocal_residual(bool request_jacobian,
218  DiffContext & context)
219 {
220  // The EigenTimeSolver always requests jacobians?
221  //libmesh_assert (request_jacobian);
222 
223  // Assemble the operator for the spatial part.
224  if (now_assembling == Matrix_A)
225  {
226  bool jacobian_computed =
227  _system.get_physics()->nonlocal_time_derivative(request_jacobian, context);
228 
229  // The user shouldn't compute a jacobian unless requested
230  libmesh_assert (request_jacobian || !jacobian_computed);
231 
232  bool jacobian_computed2 =
233  _system.get_physics()->nonlocal_constraint(jacobian_computed, context);
234 
235  // The user shouldn't compute a jacobian unless requested
236  libmesh_assert (jacobian_computed || !jacobian_computed2);
237 
238  return jacobian_computed && jacobian_computed2;
239 
240  }
241 
242  // There is now a "side" equivalent for the mass matrix
243  else if (now_assembling == Matrix_B)
244  {
245  bool mass_jacobian_computed =
246  _system.get_physics()->nonlocal_mass_residual(request_jacobian, context);
247 
248  // Scale Jacobian by -1 to get positive matrix from new negative
249  // mass_residual convention
250  context.get_elem_jacobian() *= -1.0;
251 
252  return mass_jacobian_computed;
253  }
254 
255  else
256  libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
257 }
258 
259 } // namespace libMesh
260 
261 #endif // LIBMESH_HAVE_SLEPC
bool quiet
Print extra debugging information if quiet == false.
Definition: time_solver.h:230
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on side of elem to elem_residual.
Definition: diff_physics.h:195
const DenseMatrix< Number > & get_elem_jacobian() const
Const accessor for element Jacobian.
Definition: diff_context.h:274
virtual void init() override
The initialization function.
This class provides all data required for a physics package (e.g.
Definition: diff_context.h:55
virtual bool element_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on elem to elem_residual.
Definition: diff_physics.h:125
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false) override=0
Assembles a residual in rhs and/or a jacobian in matrix, as requested.
This is a generic class that defines a solver to handle time integration of DifferentiableSystems.
Definition: time_solver.h:63
static constexpr Real TOLERANCE
virtual void solve() override
Implements the assembly of both matrices A and B, and calls the EigenSolver to compute the eigenvalue...
virtual bool nonlocal_residual(bool get_jacobian, DiffContext &) override
Forms the jacobian of the nonlocal terms.
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Adds any nonlocal time derivative contributions (e.g.
Definition: diff_physics.h:214
The libMesh namespace provides an interface to certain functionality in the library.
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Adds any nonlocal constraint contributions (e.g.
Definition: diff_physics.h:233
virtual bool mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on elem from elem_residual.
Definition: diff_physics.h:302
bool have_matrix(std::string_view mat_name) const
Definition: system.h:1888
This class provides a specific system class.
Definition: diff_system.h:54
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
T pow(const T &x)
Definition: utility.h:328
Real elem_solution_rate_derivative
The derivative of elem_solution_rate with respect to the current nonlinear solution, for use by systems with non default mass_residual terms.
Definition: diff_context.h:503
Real elem_solution_derivative
The derivative of elem_solution with respect to the current nonlinear solution.
Definition: diff_context.h:496
libmesh_assert(ctx)
unsigned int maxits
The maximum number of iterations allowed to solve the problem.
unsigned int n_iterations_reqd
After a solve, holds the number of iterations required to converge the requested number of eigenpairs...
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:181
virtual ~EigenTimeSolver()
Destructor.
NowAssembling now_assembling
Flag which controls the internals of element_residual() and side_residual().
double tol
The linear solver tolerance to be used when solving the eigenvalue problem.
OStreamProxy out
SparseMatrix< Number > * matrix
The system matrix.
unsigned int n_basis_vectors_to_use
The number of basis vectors to use in the computation.
EigenTimeSolver(sys_type &s)
Constructor.
virtual bool side_residual(bool get_jacobian, DiffContext &) override
Forms the jacobian of the boundary terms.
virtual void reinit() override
The reinitialization function.
virtual bool side_mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on side of elem from elem_residual.
Definition: diff_physics.h:320
virtual bool side_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on side of elem to elem_residual.
Definition: diff_physics.h:174
virtual bool element_residual(bool get_jacobian, DiffContext &) override
Forms either the spatial (Jacobian) or mass matrix part of the operator, depending on which is reques...
unsigned int n_converged_eigenpairs
After a solve, holds the number of eigenpairs successfully converged.
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on elem to elem_residual.
Definition: diff_physics.h:144
SparseMatrix< Number > & add_matrix(std::string_view mat_name, ParallelType type=PARALLEL, MatrixBuildType mat_build_type=MatrixBuildType::AUTOMATIC)
Adds the additional matrix mat_name to this system.
Definition: system.C:997
const SparseMatrix< Number > & get_matrix(std::string_view mat_name) const
Definition: system.C:1102
unsigned int n_eigenpairs_to_compute
The number of eigenvectors/values to be computed.
The matrix associated with the spatial part of the operator.
const SparseMatrix< Number > & get_system_matrix() const
The matrix associated with the time derivative (mass matrix).
virtual bool nonlocal_mass_residual(bool request_jacobian, DiffContext &c)
Subtracts any nonlocal mass vector contributions (e.g.
Definition: diff_physics.C:57
std::unique_ptr< EigenSolver< Number > > eigen_solver
The EigenSolver object.
This class provides an interface to solvers for eigenvalue problems.
Definition: eigen_solver.h:58