libMesh
Public Types | Public Member Functions | Static Public Member Functions | Public Attributes | Protected Types | Protected Member Functions | Protected Attributes | Static Protected Attributes | Private Types | Private Attributes | List of all members
libMesh::EigenTimeSolver Class Reference

The name of this class is confusing...it's meant to refer to the base class (TimeSolver) while still telling one that it's for solving (generalized) EigenValue problems that arise from finite element discretizations. More...

#include <eigen_time_solver.h>

Inheritance diagram for libMesh::EigenTimeSolver:
[legend]

Public Types

typedef DifferentiableSystem sys_type
 The type of system. More...
 
typedef TimeSolver Parent
 The parent class. More...
 

Public Member Functions

 EigenTimeSolver (sys_type &s)
 Constructor. More...
 
virtual ~EigenTimeSolver ()
 Destructor. More...
 
virtual void init () override
 The initialization function. More...
 
virtual void reinit () override
 The reinitialization function. More...
 
virtual void solve () override
 Implements the assembly of both matrices A and B, and calls the EigenSolver to compute the eigenvalues. More...
 
virtual void advance_timestep () override
 It doesn't make sense to advance the timestep, so we shouldn't call this. More...
 
Real error_order () const
 error convergence order against deltat is not applicable to an eigenvalue problem. More...
 
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 requested. More...
 
virtual bool side_residual (bool get_jacobian, DiffContext &) override
 Forms the jacobian of the boundary terms. More...
 
virtual bool nonlocal_residual (bool get_jacobian, DiffContext &) override
 Forms the jacobian of the nonlocal terms. More...
 
virtual Real du (const SystemNorm &) const override
 
virtual bool is_steady () const override
 This is effectively a steady-state solver. More...
 
virtual void init_data ()
 The data initialization function. More...
 
virtual void adjoint_advance_timestep ()
 This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. More...
 
virtual void retrieve_timestep ()
 This method retrieves all the stored solutions at the current system.time. More...
 
virtual void before_timestep ()
 This method is for subclasses or users to override to do arbitrary processing between timesteps. More...
 
const sys_typesystem () const
 
sys_typesystem ()
 
virtual std::unique_ptr< DiffSolver > & diff_solver ()
 An implicit linear or nonlinear solver to use at each timestep. More...
 
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver ()
 An implicit linear solver to use for adjoint and sensitivity problems. More...
 
void set_solution_history (const SolutionHistory &_solution_history)
 A setter function users will employ if they need to do something other than save no solution history. More...
 
bool is_adjoint () const
 Accessor for querying whether we need to do a primal or adjoint solve. More...
 
void set_is_adjoint (bool _is_adjoint_value)
 Accessor for setting whether we need to do a primal or adjoint solve. More...
 

Static Public Member Functions

static std::string get_info ()
 Gets a string containing the reference information. More...
 
static void print_info (std::ostream &out=libMesh::out)
 Prints the reference information, by default to libMesh::out. More...
 
static unsigned int n_objects ()
 Prints the number of outstanding (created, but not yet destroyed) objects. More...
 
static void enable_print_counter_info ()
 Methods to enable/disable the reference counter output from print_info() More...
 
static void disable_print_counter_info ()
 

Public Attributes

std::unique_ptr< EigenSolver< Number > > eigen_solver
 The EigenSolver object. More...
 
Real tol
 The linear solver tolerance to be used when solving the eigenvalue problem. More...
 
unsigned int maxits
 The maximum number of iterations allowed to solve the problem. More...
 
unsigned int n_eigenpairs_to_compute
 The number of eigenvectors/values to be computed. More...
 
unsigned int n_basis_vectors_to_use
 The number of basis vectors to use in the computation. More...
 
unsigned int n_converged_eigenpairs
 After a solve, holds the number of eigenpairs successfully converged. More...
 
unsigned int n_iterations_reqd
 After a solve, holds the number of iterations required to converge the requested number of eigenpairs. More...
 
bool quiet
 Print extra debugging information if quiet == false. More...
 
unsigned int reduce_deltat_on_diffsolver_failure
 This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat and let the DiffSolver repeat the latest failed solve with a reduced timestep. More...
 

Protected Types

typedef bool(DifferentiablePhysics::* ResFuncType) (bool, DiffContext &)
 Definitions of argument types for use in refactoring subclasses. More...
 
typedef void(DiffContext::* ReinitFuncType) (Real)
 
typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 Data structure to log the information. More...
 

Protected Member Functions

void increment_constructor_count (const std::string &name)
 Increments the construction counter. More...
 
void increment_destructor_count (const std::string &name)
 Increments the destruction counter. More...
 

Protected Attributes

std::unique_ptr< DiffSolver_diff_solver
 An implicit linear or nonlinear solver to use at each timestep. More...
 
std::unique_ptr< LinearSolver< Number > > _linear_solver
 An implicit linear solver to use for adjoint problems. More...
 
sys_type_system
 A reference to the system we are solving. More...
 
std::unique_ptr< SolutionHistorysolution_history
 A std::unique_ptr to a SolutionHistory object. More...
 

Static Protected Attributes

static Counts _counts
 Actually holds the data. More...
 
static Threads::atomic< unsigned int_n_objects
 The number of objects. More...
 
static Threads::spin_mutex _mutex
 Mutual exclusion object to enable thread-safe reference counting. More...
 
static bool _enable_print_counter = true
 Flag to control whether reference count information is printed when print_info is called. More...
 

Private Types

enum  NowAssembling { Matrix_A, Matrix_B, Invalid_Matrix }
 

Private Attributes

NowAssembling now_assembling
 Flag which controls the internals of element_residual() and side_residual(). More...
 
bool _is_adjoint
 This boolean tells the TimeSolver whether we are solving a primal or adjoint problem. More...
 

Detailed Description

The name of this class is confusing...it's meant to refer to the base class (TimeSolver) while still telling one that it's for solving (generalized) EigenValue problems that arise from finite element discretizations.

For a time-dependent problem du/dt=F(u), with a steady solution 0=F(u_0), we look at the time evolution of a small perturbation, p=u-u_0, for which the (linearized) governing equation is

dp/dt = F'(u_0)p

where F'(u_0) is the Jacobian. The generalized eigenvalue problem arises by considering perturbations of the general form p = exp(lambda*t)x, which leads to

Ax = lambda*Bx

where A is the (discretized by FEM) Jacobian matrix and B is the (discretized by FEM) mass matrix.

The EigenSystem class (by Steffen Petersen) is related but does not fall under the FEMSystem paradigm invented by Roy Stogner. The EigenSolver class (also by Steffen) is meant to provide a generic "linear solver" interface for EigenValue software. The only current concrete implementation is a SLEPc-based eigensolver class, which we make use of here as well.

Author
John W. Peterson
Date
2007

Definition at line 66 of file eigen_time_solver.h.

Member Typedef Documentation

◆ Counts

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts
protectedinherited

Data structure to log the information.

The log is identified by the class name.

Definition at line 117 of file reference_counter.h.

◆ Parent

The parent class.

Definition at line 77 of file eigen_time_solver.h.

◆ ReinitFuncType

typedef void(DiffContext::* libMesh::TimeSolver::ReinitFuncType) (Real)
protectedinherited

Definition at line 272 of file time_solver.h.

◆ ResFuncType

typedef bool(DifferentiablePhysics::* libMesh::TimeSolver::ResFuncType) (bool, DiffContext &)
protectedinherited

Definitions of argument types for use in refactoring subclasses.

Definition at line 270 of file time_solver.h.

◆ sys_type

The type of system.

Definition at line 72 of file eigen_time_solver.h.

Member Enumeration Documentation

◆ NowAssembling

Enumerator
Matrix_A 

The matrix associated with the spatial part of the operator.

Matrix_B 

The matrix associated with the time derivative (mass matrix).

Invalid_Matrix 

The enum is in an invalid state.

Definition at line 198 of file eigen_time_solver.h.

198  {
202  Matrix_A,
203 
207  Matrix_B,
208 
213  };

Constructor & Destructor Documentation

◆ EigenTimeSolver()

libMesh::EigenTimeSolver::EigenTimeSolver ( sys_type s)
explicit

Constructor.

Requires a reference to the system to be solved.

Definition at line 33 of file eigen_time_solver.C.

34  : Parent(s),
36  tol(pow(TOLERANCE, 5./3.)),
37  maxits(1000),
42 {
43  libmesh_experimental();
44  eigen_solver->set_eigenproblem_type(GHEP);//or GNHEP
45  eigen_solver->set_position_of_spectrum(LARGEST_MAGNITUDE);
46 }

References eigen_solver, libMesh::GHEP, and libMesh::LARGEST_MAGNITUDE.

◆ ~EigenTimeSolver()

libMesh::EigenTimeSolver::~EigenTimeSolver ( )
virtual

Destructor.

Definition at line 48 of file eigen_time_solver.C.

49 {
50 }

Member Function Documentation

◆ adjoint_advance_timestep()

void libMesh::TimeSolver::adjoint_advance_timestep ( )
virtualinherited

This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed.

This will be done before every UnsteadySolver::adjoint_solve().

Reimplemented in libMesh::UnsteadySolver, and libMesh::NewmarkSolver.

Definition at line 105 of file time_solver.C.

106 {
107 }

◆ advance_timestep()

virtual void libMesh::EigenTimeSolver::advance_timestep ( )
inlineoverridevirtual

It doesn't make sense to advance the timestep, so we shouldn't call this.

Reimplemented from libMesh::TimeSolver.

Definition at line 112 of file eigen_time_solver.h.

112 {}

◆ before_timestep()

virtual void libMesh::TimeSolver::before_timestep ( )
inlinevirtualinherited

This method is for subclasses or users to override to do arbitrary processing between timesteps.

Definition at line 166 of file time_solver.h.

166 {}

◆ diff_solver()

virtual std::unique_ptr<DiffSolver>& libMesh::TimeSolver::diff_solver ( )
inlinevirtualinherited

An implicit linear or nonlinear solver to use at each timestep.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 181 of file time_solver.h.

181 { return _diff_solver; }

References libMesh::TimeSolver::_diff_solver.

Referenced by libMesh::TimeSolver::init(), libMesh::TimeSolver::init_data(), libMesh::TimeSolver::reinit(), and libMesh::TimeSolver::solve().

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

107 {
108  _enable_print_counter = false;
109  return;
110 }

References libMesh::ReferenceCounter::_enable_print_counter.

Referenced by libMesh::LibMeshInit::LibMeshInit().

◆ du()

virtual Real libMesh::EigenTimeSolver::du ( const SystemNorm ) const
inlineoverridevirtual
Returns
0, but derived classes should override this function to compute the size of the difference between successive solution iterates ||u^{n+1} - u^{n}|| in some norm.

Implements libMesh::TimeSolver.

Definition at line 144 of file eigen_time_solver.h.

144 { return 0.; }

◆ element_residual()

bool libMesh::EigenTimeSolver::element_residual ( bool  get_jacobian,
DiffContext context 
)
overridevirtual

Forms either the spatial (Jacobian) or mass matrix part of the operator, depending on which is requested.

Implements libMesh::TimeSolver.

Definition at line 128 of file eigen_time_solver.C.

130 {
131  // The EigenTimeSolver always computes jacobians!
132  libmesh_assert (request_jacobian);
133 
134  context.elem_solution_rate_derivative = 1;
135  context.elem_solution_derivative = 1;
136 
137  // Assemble the operator for the spatial part.
138  if (now_assembling == Matrix_A)
139  {
140  bool jacobian_computed =
141  _system.get_physics()->element_time_derivative(request_jacobian, context);
142 
143  // The user shouldn't compute a jacobian unless requested
144  libmesh_assert(request_jacobian || !jacobian_computed);
145 
146  bool jacobian_computed2 =
147  _system.get_physics()->element_constraint(jacobian_computed, context);
148 
149  // The user shouldn't compute a jacobian unless requested
150  libmesh_assert (jacobian_computed || !jacobian_computed2);
151 
152  return jacobian_computed && jacobian_computed2;
153 
154  }
155 
156  // Assemble the mass matrix operator
157  else if (now_assembling == Matrix_B)
158  {
159  bool mass_jacobian_computed =
160  _system.get_physics()->mass_residual(request_jacobian, context);
161 
162  // Scale Jacobian by -1 to get positive matrix from new negative
163  // mass_residual convention
164  context.get_elem_jacobian() *= -1.0;
165 
166  return mass_jacobian_computed;
167  }
168 
169  else
170  libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
171 }

References libMesh::TimeSolver::_system, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DifferentiablePhysics::element_constraint(), libMesh::DifferentiablePhysics::element_time_derivative(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DifferentiableSystem::get_physics(), libMesh::libmesh_assert(), libMesh::DifferentiablePhysics::mass_residual(), Matrix_A, Matrix_B, and now_assembling.

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

101 {
102  _enable_print_counter = true;
103  return;
104 }

References libMesh::ReferenceCounter::_enable_print_counter.

◆ error_order()

Real libMesh::EigenTimeSolver::error_order ( ) const
inline

error convergence order against deltat is not applicable to an eigenvalue problem.

Definition at line 118 of file eigen_time_solver.h.

118 { return 0.; }

◆ get_info()

std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (const auto & pr : _counts)
59  {
60  const std::string name(pr.first);
61  const unsigned int creations = pr.second.first;
62  const unsigned int destructions = pr.second.second;
63 
64  oss << "| " << name << " reference count information:\n"
65  << "| Creations: " << creations << '\n'
66  << "| Destructions: " << destructions << '\n';
67  }
68 
69  oss << " ---------------------------------------------------------------------------- \n";
70 
71  return oss.str();
72 
73 #else
74 
75  return "";
76 
77 #endif
78 }

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

◆ increment_constructor_count()

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the construction counter.

Should be called in the constructor of any derived class that will be reference counted.

Definition at line 181 of file reference_counter.h.

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

◆ increment_destructor_count()

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

Increments the destruction counter.

Should be called in the destructor of any derived class that will be reference counted.

Definition at line 194 of file reference_counter.h.

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }

References libMesh::ReferenceCounter::_counts, libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

◆ init()

void libMesh::EigenTimeSolver::init ( )
overridevirtual

The initialization function.

This method is used to initialize internal data structures before a simulation begins.

Reimplemented from libMesh::TimeSolver.

Definition at line 57 of file eigen_time_solver.C.

58 {
59  // Add matrix "B" to _system if not already there.
60  // The user may have already added a matrix "B" before
61  // calling the System initialization. This would be
62  // necessary if e.g. the System originally started life
63  // with a different type of TimeSolver and only later
64  // had its TimeSolver changed to an EigenTimeSolver.
65  if (!_system.have_matrix("B"))
66  _system.add_matrix("B");
67 }

References libMesh::TimeSolver::_system, libMesh::ImplicitSystem::add_matrix(), and libMesh::ImplicitSystem::have_matrix().

◆ init_data()

void libMesh::TimeSolver::init_data ( )
virtualinherited

The data initialization function.

This method is used to initialize internal data structures after the underlying System has been initialized

Reimplemented in libMesh::UnsteadySolver, and libMesh::SecondOrderUnsteadySolver.

Definition at line 76 of file time_solver.C.

77 {
78  this->diff_solver()->init();
79 
80  if (libMesh::on_command_line("--solver-system-names"))
81  this->linear_solver()->init((_system.name()+"_").c_str());
82  else
83  this->linear_solver()->init();
84 }

References libMesh::TimeSolver::_system, libMesh::TimeSolver::diff_solver(), libMesh::TimeSolver::linear_solver(), libMesh::System::name(), and libMesh::on_command_line().

Referenced by libMesh::UnsteadySolver::init_data().

◆ is_adjoint()

bool libMesh::TimeSolver::is_adjoint ( ) const
inlineinherited

Accessor for querying whether we need to do a primal or adjoint solve.

Definition at line 232 of file time_solver.h.

233  { return _is_adjoint; }

References libMesh::TimeSolver::_is_adjoint.

Referenced by libMesh::FEMSystem::build_context().

◆ is_steady()

virtual bool libMesh::EigenTimeSolver::is_steady ( ) const
inlineoverridevirtual

This is effectively a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 149 of file eigen_time_solver.h.

149 { return true; }

◆ linear_solver()

virtual std::unique_ptr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver ( )
inlinevirtualinherited

An implicit linear solver to use for adjoint and sensitivity problems.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 186 of file time_solver.h.

186 { return _linear_solver; }

References libMesh::TimeSolver::_linear_solver.

Referenced by libMesh::TimeSolver::init(), libMesh::TimeSolver::init_data(), and libMesh::TimeSolver::reinit().

◆ n_objects()

static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 83 of file reference_counter.h.

84  { return _n_objects; }

References libMesh::ReferenceCounter::_n_objects.

◆ nonlocal_residual()

bool libMesh::EigenTimeSolver::nonlocal_residual ( bool  get_jacobian,
DiffContext context 
)
overridevirtual

Forms the jacobian of the nonlocal terms.

Implements libMesh::TimeSolver.

Definition at line 222 of file eigen_time_solver.C.

224 {
225  // The EigenTimeSolver always requests jacobians?
226  //libmesh_assert (request_jacobian);
227 
228  // Assemble the operator for the spatial part.
229  if (now_assembling == Matrix_A)
230  {
231  bool jacobian_computed =
232  _system.get_physics()->nonlocal_time_derivative(request_jacobian, context);
233 
234  // The user shouldn't compute a jacobian unless requested
235  libmesh_assert (request_jacobian || !jacobian_computed);
236 
237  bool jacobian_computed2 =
238  _system.get_physics()->nonlocal_constraint(jacobian_computed, context);
239 
240  // The user shouldn't compute a jacobian unless requested
241  libmesh_assert (jacobian_computed || !jacobian_computed2);
242 
243  return jacobian_computed && jacobian_computed2;
244 
245  }
246 
247  // There is now a "side" equivalent for the mass matrix
248  else if (now_assembling == Matrix_B)
249  {
250  bool mass_jacobian_computed =
251  _system.get_physics()->nonlocal_mass_residual(request_jacobian, context);
252 
253  // Scale Jacobian by -1 to get positive matrix from new negative
254  // mass_residual convention
255  context.get_elem_jacobian() *= -1.0;
256 
257  return mass_jacobian_computed;
258  }
259 
260  else
261  libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
262 }

References libMesh::TimeSolver::_system, libMesh::DiffContext::get_elem_jacobian(), libMesh::DifferentiableSystem::get_physics(), libMesh::libmesh_assert(), Matrix_A, Matrix_B, libMesh::DifferentiablePhysics::nonlocal_constraint(), libMesh::DifferentiablePhysics::nonlocal_mass_residual(), libMesh::DifferentiablePhysics::nonlocal_time_derivative(), and now_assembling.

◆ print_info()

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

Prints the reference information, by default to libMesh::out.

Definition at line 87 of file reference_counter.C.

88 {
90  out_stream << ReferenceCounter::get_info();
91 }

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

◆ reinit()

void libMesh::EigenTimeSolver::reinit ( )
overridevirtual

The reinitialization function.

This method is used after changes in the mesh

Reimplemented from libMesh::TimeSolver.

Definition at line 52 of file eigen_time_solver.C.

53 {
54  // empty...
55 }

◆ retrieve_timestep()

void libMesh::TimeSolver::retrieve_timestep ( )
virtualinherited

This method retrieves all the stored solutions at the current system.time.

Reimplemented in libMesh::UnsteadySolver, and libMesh::SecondOrderUnsteadySolver.

Definition at line 109 of file time_solver.C.

110 {
111 }

◆ set_is_adjoint()

void libMesh::TimeSolver::set_is_adjoint ( bool  _is_adjoint_value)
inlineinherited

Accessor for setting whether we need to do a primal or adjoint solve.

Definition at line 239 of file time_solver.h.

240  { _is_adjoint = _is_adjoint_value; }

References libMesh::TimeSolver::_is_adjoint.

Referenced by libMesh::DifferentiableSystem::adjoint_solve(), libMesh::FEMSystem::postprocess(), and libMesh::DifferentiableSystem::solve().

◆ set_solution_history()

void libMesh::TimeSolver::set_solution_history ( const SolutionHistory _solution_history)
inherited

A setter function users will employ if they need to do something other than save no solution history.

Definition at line 96 of file time_solver.C.

97 {
98  solution_history = _solution_history.clone();
99 }

References libMesh::SolutionHistory::clone(), and libMesh::TimeSolver::solution_history.

◆ side_residual()

bool libMesh::EigenTimeSolver::side_residual ( bool  get_jacobian,
DiffContext context 
)
overridevirtual

Forms the jacobian of the boundary terms.

Implements libMesh::TimeSolver.

Definition at line 175 of file eigen_time_solver.C.

177 {
178  // The EigenTimeSolver always requests jacobians?
179  //libmesh_assert (request_jacobian);
180 
181  context.elem_solution_rate_derivative = 1;
182  context.elem_solution_derivative = 1;
183 
184  // Assemble the operator for the spatial part.
185  if (now_assembling == Matrix_A)
186  {
187  bool jacobian_computed =
188  _system.get_physics()->side_time_derivative(request_jacobian, context);
189 
190  // The user shouldn't compute a jacobian unless requested
191  libmesh_assert (request_jacobian || !jacobian_computed);
192 
193  bool jacobian_computed2 =
194  _system.get_physics()->side_constraint(jacobian_computed, context);
195 
196  // The user shouldn't compute a jacobian unless requested
197  libmesh_assert (jacobian_computed || !jacobian_computed2);
198 
199  return jacobian_computed && jacobian_computed2;
200 
201  }
202 
203  // There is now a "side" equivalent for the mass matrix
204  else if (now_assembling == Matrix_B)
205  {
206  bool mass_jacobian_computed =
207  _system.get_physics()->side_mass_residual(request_jacobian, context);
208 
209  // Scale Jacobian by -1 to get positive matrix from new negative
210  // mass_residual convention
211  context.get_elem_jacobian() *= -1.0;
212 
213  return mass_jacobian_computed;
214  }
215 
216  else
217  libmesh_error_msg("Unrecognized value now_assembling = " << now_assembling);
218 }

References libMesh::TimeSolver::_system, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DiffContext::get_elem_jacobian(), libMesh::DifferentiableSystem::get_physics(), libMesh::libmesh_assert(), Matrix_A, Matrix_B, now_assembling, libMesh::DifferentiablePhysics::side_constraint(), libMesh::DifferentiablePhysics::side_mass_residual(), and libMesh::DifferentiablePhysics::side_time_derivative().

◆ solve()

void libMesh::EigenTimeSolver::solve ( )
overridevirtual

Implements the assembly of both matrices A and B, and calls the EigenSolver to compute the eigenvalues.

Reimplemented from libMesh::TimeSolver.

Definition at line 69 of file eigen_time_solver.C.

70 {
71  // The standard implementation is basically to call:
72  // _diff_solver->solve();
73  // which internally assembles (when necessary) matrices and vectors
74  // and calls linear solver software while also doing Newton steps (see newton_solver.C)
75  //
76  // The element_residual and side_residual functions below control
77  // what happens in the interior of the element assembly loops.
78  // We have a system reference, so it's possible to call _system.assembly()
79  // ourselves if we want to...
80  //
81  // Interestingly, for the EigenSolver we don't need residuals...just Jacobians.
82  // The Jacobian should therefore always be requested, and always return
83  // jacobian_computed as being true.
84 
85  // The basic plan of attack is:
86  // .) Construct the Jacobian using _system.assembly(true,true) as we
87  // would for a steady system. Use a flag in this class to
88  // control behavior in element_residual and side_residual
89  // .) Swap _system.matrix to matrix "B" (be sure to add this extra matrix during init)
90  // .) Call _system.assembly(true,true) again, using the flag in element_residual
91  // and side_residual to only get the mass matrix terms.
92  // .) Send A and B to Steffen's EigenSolver interface.
93 
94  // Assemble the spatial part (matrix A) of the operator
95  if (!this->quiet)
96  libMesh::out << "Assembling matrix A." << std::endl;
97  _system.matrix = &( _system.get_matrix ("System Matrix") );
98  this->now_assembling = Matrix_A;
99  _system.assembly(true, true);
100  //_system.matrix->print_matlab("matrix_A.m");
101 
102  // Point the system's matrix at B, call assembly again.
103  if (!this->quiet)
104  libMesh::out << "Assembling matrix B." << std::endl;
105  _system.matrix = &( _system.get_matrix ("B") );
106  this->now_assembling = Matrix_B;
107  _system.assembly(true, true);
108  //_system.matrix->print_matlab("matrix_B.m");
109 
110  // Send matrices A, B to Steffen's SlepcEigenSolver interface
111  //libmesh_here();
112  if (!this->quiet)
113  libMesh::out << "Calling the EigenSolver." << std::endl;
114  std::pair<unsigned int, unsigned int> solve_data =
115  eigen_solver->solve_generalized (_system.get_matrix ("System Matrix"),
116  _system.get_matrix ("B"),
119  tol,
120  maxits);
121 
122  this->n_converged_eigenpairs = solve_data.first;
123  this->n_iterations_reqd = solve_data.second;
124 }

References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::assembly(), eigen_solver, libMesh::ImplicitSystem::get_matrix(), libMesh::ImplicitSystem::matrix, Matrix_A, Matrix_B, maxits, n_basis_vectors_to_use, n_converged_eigenpairs, n_eigenpairs_to_compute, n_iterations_reqd, now_assembling, libMesh::out, libMesh::TimeSolver::quiet, and tol.

◆ system() [1/2]

sys_type& libMesh::TimeSolver::system ( )
inlineinherited
Returns
A writable reference to the system we are solving.

Definition at line 176 of file time_solver.h.

176 { return _system; }

References libMesh::TimeSolver::_system.

◆ system() [2/2]

const sys_type& libMesh::TimeSolver::system ( ) const
inlineinherited
Returns
A constant reference to the system we are solving.

Definition at line 171 of file time_solver.h.

171 { return _system; }

References libMesh::TimeSolver::_system.

Referenced by libMesh::TimeSolver::reinit(), and libMesh::TimeSolver::solve().

Member Data Documentation

◆ _counts

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited

◆ _diff_solver

std::unique_ptr<DiffSolver> libMesh::TimeSolver::_diff_solver
protectedinherited

An implicit linear or nonlinear solver to use at each timestep.

Definition at line 247 of file time_solver.h.

Referenced by libMesh::NewmarkSolver::compute_initial_accel(), libMesh::TimeSolver::diff_solver(), and libMesh::UnsteadySolver::solve().

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

Flag to control whether reference count information is printed when print_info is called.

Definition at line 141 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

◆ _is_adjoint

bool libMesh::TimeSolver::_is_adjoint
privateinherited

This boolean tells the TimeSolver whether we are solving a primal or adjoint problem.

Definition at line 280 of file time_solver.h.

Referenced by libMesh::TimeSolver::is_adjoint(), and libMesh::TimeSolver::set_is_adjoint().

◆ _linear_solver

std::unique_ptr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver
protectedinherited

An implicit linear solver to use for adjoint problems.

Definition at line 252 of file time_solver.h.

Referenced by libMesh::TimeSolver::linear_solver().

◆ _mutex

Threads::spin_mutex libMesh::ReferenceCounter::_mutex
staticprotectedinherited

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects
staticprotectedinherited

The number of objects.

Print the reference count information when the number returns to 0.

Definition at line 130 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

◆ _system

sys_type& libMesh::TimeSolver::_system
protectedinherited

A reference to the system we are solving.

Definition at line 257 of file time_solver.h.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::NewmarkSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::NewmarkSolver::compute_initial_accel(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::UnsteadySolver::du(), libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), element_residual(), libMesh::SecondOrderUnsteadySolver::init(), libMesh::UnsteadySolver::init(), libMesh::TimeSolver::init(), init(), libMesh::SecondOrderUnsteadySolver::init_data(), libMesh::UnsteadySolver::init_data(), libMesh::TimeSolver::init_data(), libMesh::EulerSolver::nonlocal_residual(), libMesh::Euler2Solver::nonlocal_residual(), nonlocal_residual(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::SecondOrderUnsteadySolver::old_solution_accel(), libMesh::SecondOrderUnsteadySolver::old_solution_rate(), libMesh::NewmarkSolver::project_initial_accel(), libMesh::SecondOrderUnsteadySolver::project_initial_rate(), libMesh::SecondOrderUnsteadySolver::reinit(), libMesh::UnsteadySolver::reinit(), libMesh::TimeSolver::reinit(), libMesh::UnsteadySolver::retrieve_timestep(), side_residual(), libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), solve(), and libMesh::TimeSolver::system().

◆ eigen_solver

std::unique_ptr<EigenSolver<Number> > libMesh::EigenTimeSolver::eigen_solver

The EigenSolver object.

This is what actually makes the calls to SLEPc.

Definition at line 155 of file eigen_time_solver.h.

Referenced by EigenTimeSolver(), and solve().

◆ maxits

unsigned int libMesh::EigenTimeSolver::maxits

The maximum number of iterations allowed to solve the problem.

Definition at line 166 of file eigen_time_solver.h.

Referenced by solve().

◆ n_basis_vectors_to_use

unsigned int libMesh::EigenTimeSolver::n_basis_vectors_to_use

The number of basis vectors to use in the computation.

According to ex16, the number of basis vectors must be >= the number of eigenpairs requested, and ncv >= 2*nev is recommended. Increasing this number, even by a little bit, can greatly reduce the number of (EigenSolver) iterations required to compute the desired number of eigenpairs, but the cost per iteration goes up drastically as well.

Definition at line 182 of file eigen_time_solver.h.

Referenced by solve().

◆ n_converged_eigenpairs

unsigned int libMesh::EigenTimeSolver::n_converged_eigenpairs

After a solve, holds the number of eigenpairs successfully converged.

Definition at line 188 of file eigen_time_solver.h.

Referenced by solve().

◆ n_eigenpairs_to_compute

unsigned int libMesh::EigenTimeSolver::n_eigenpairs_to_compute

The number of eigenvectors/values to be computed.

Definition at line 171 of file eigen_time_solver.h.

Referenced by solve().

◆ n_iterations_reqd

unsigned int libMesh::EigenTimeSolver::n_iterations_reqd

After a solve, holds the number of iterations required to converge the requested number of eigenpairs.

Definition at line 194 of file eigen_time_solver.h.

Referenced by solve().

◆ now_assembling

NowAssembling libMesh::EigenTimeSolver::now_assembling
private

Flag which controls the internals of element_residual() and side_residual().

Definition at line 218 of file eigen_time_solver.h.

Referenced by element_residual(), nonlocal_residual(), side_residual(), and solve().

◆ quiet

bool libMesh::TimeSolver::quiet
inherited

Print extra debugging information if quiet == false.

Definition at line 191 of file time_solver.h.

Referenced by libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), and solve().

◆ reduce_deltat_on_diffsolver_failure

unsigned int libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure
inherited

This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat and let the DiffSolver repeat the latest failed solve with a reduced timestep.

Note
This has no effect for SteadySolvers.
You must set at least one of the DiffSolver flags "continue_after_max_iterations" or "continue_after_backtrack_failure" to allow the TimeSolver to retry the solve.

Definition at line 220 of file time_solver.h.

Referenced by libMesh::TwostepTimeSolver::solve(), and libMesh::UnsteadySolver::solve().

◆ solution_history

std::unique_ptr<SolutionHistory> libMesh::TimeSolver::solution_history
protectedinherited

A std::unique_ptr to a SolutionHistory object.

Default is NoSolutionHistory, which the user can override by declaring a different kind of SolutionHistory in the application

Definition at line 264 of file time_solver.h.

Referenced by libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::UnsteadySolver::retrieve_timestep(), and libMesh::TimeSolver::set_solution_history().

◆ tol

Real libMesh::EigenTimeSolver::tol

The linear solver tolerance to be used when solving the eigenvalue problem.

FIXME: need more info...

Definition at line 161 of file eigen_time_solver.h.

Referenced by solve().


The documentation for this class was generated from the following files:
libMesh::DifferentiablePhysics::nonlocal_mass_residual
virtual bool nonlocal_mass_residual(bool request_jacobian, DiffContext &c)
Subtracts any nonlocal mass vector contributions (e.g.
Definition: diff_physics.C:63
libMesh::DifferentiablePhysics::side_constraint
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on side of elem to elem_residual.
Definition: diff_physics.h:192
libMesh::TimeSolver::quiet
bool quiet
Print extra debugging information if quiet == false.
Definition: time_solver.h:191
libMesh::TimeSolver::_diff_solver
std::unique_ptr< DiffSolver > _diff_solver
An implicit linear or nonlinear solver to use at each timestep.
Definition: time_solver.h:247
libMesh::EigenTimeSolver::n_basis_vectors_to_use
unsigned int n_basis_vectors_to_use
The number of basis vectors to use in the computation.
Definition: eigen_time_solver.h:182
libMesh::EigenTimeSolver::tol
Real tol
The linear solver tolerance to be used when solving the eigenvalue problem.
Definition: eigen_time_solver.h:161
libMesh::ImplicitSystem::add_matrix
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
Adds the additional matrix mat_name to this system.
Definition: implicit_system.C:202
libMesh::DifferentiablePhysics::nonlocal_constraint
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Adds any nonlocal constraint contributions (e.g.
Definition: diff_physics.h:228
libMesh::ReferenceCounter::_counts
static Counts _counts
Actually holds the data.
Definition: reference_counter.h:122
libMesh::TOLERANCE
static const Real TOLERANCE
Definition: libmesh_common.h:128
libMesh::EigenTimeSolver::maxits
unsigned int maxits
The maximum number of iterations allowed to solve the problem.
Definition: eigen_time_solver.h:166
libMesh::ReferenceCounter::_n_objects
static Threads::atomic< unsigned int > _n_objects
The number of objects.
Definition: reference_counter.h:130
libMesh::EigenTimeSolver::Parent
TimeSolver Parent
The parent class.
Definition: eigen_time_solver.h:77
libMesh::ReferenceCounter::get_info
static std::string get_info()
Gets a string containing the reference information.
Definition: reference_counter.C:47
libMesh::DifferentiablePhysics::nonlocal_time_derivative
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Adds any nonlocal time derivative contributions (e.g.
Definition: diff_physics.h:210
libMesh::TimeSolver::_is_adjoint
bool _is_adjoint
This boolean tells the TimeSolver whether we are solving a primal or adjoint problem.
Definition: time_solver.h:280
libMesh::DifferentiablePhysics::element_constraint
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on elem to elem_residual.
Definition: diff_physics.h:143
libMesh::libmesh_assert
libmesh_assert(ctx)
libMesh::TimeSolver::linear_solver
virtual std::unique_ptr< LinearSolver< Number > > & linear_solver()
An implicit linear solver to use for adjoint and sensitivity problems.
Definition: time_solver.h:186
libMesh::EigenTimeSolver::Matrix_B
The matrix associated with the time derivative (mass matrix).
Definition: eigen_time_solver.h:207
libMesh::DifferentiablePhysics::element_time_derivative
virtual bool element_time_derivative(bool request_jacobian, DiffContext &)
Adds the time derivative contribution on elem to elem_residual.
Definition: diff_physics.h:125
libMesh::Threads::spin_mtx
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:29
libMesh::EigenTimeSolver::n_converged_eigenpairs
unsigned int n_converged_eigenpairs
After a solve, holds the number of eigenpairs successfully converged.
Definition: eigen_time_solver.h:188
libMesh::TimeSolver::diff_solver
virtual std::unique_ptr< DiffSolver > & diff_solver()
An implicit linear or nonlinear solver to use at each timestep.
Definition: time_solver.h:181
libMesh::ImplicitSystem::matrix
SparseMatrix< Number > * matrix
The system matrix.
Definition: implicit_system.h:393
libMesh::DifferentiablePhysics::side_mass_residual
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:336
libMesh::EigenTimeSolver::n_iterations_reqd
unsigned int n_iterations_reqd
After a solve, holds the number of iterations required to converge the requested number of eigenpairs...
Definition: eigen_time_solver.h:194
libMesh::EigenTimeSolver::Invalid_Matrix
The enum is in an invalid state.
Definition: eigen_time_solver.h:212
libMesh::ImplicitSystem::have_matrix
bool have_matrix(const std::string &mat_name) const
Definition: implicit_system.h:439
std::pow
double pow(double a, int b)
Definition: libmesh_augment_std_namespace.h:58
libMesh::TimeSolver::_system
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:257
libMesh::LARGEST_MAGNITUDE
Definition: enum_eigen_solver_type.h:75
libMesh::System::name
const std::string & name() const
Definition: system.h:2067
libMesh::EigenTimeSolver::now_assembling
NowAssembling now_assembling
Flag which controls the internals of element_residual() and side_residual().
Definition: eigen_time_solver.h:218
libMesh::GHEP
Definition: enum_eigen_solver_type.h:58
libMesh::DifferentiableSystem::get_physics
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:181
libMesh::on_command_line
bool on_command_line(std::string arg)
Definition: libmesh.C:898
libMesh::DifferentiableSystem::assembly
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.
libMesh::EigenTimeSolver::n_eigenpairs_to_compute
unsigned int n_eigenpairs_to_compute
The number of eigenvectors/values to be computed.
Definition: eigen_time_solver.h:171
libMesh::EigenSolver::build
static std::unique_ptr< EigenSolver< T > > build(const Parallel::Communicator &comm_in, const SolverPackage solver_package=SLEPC_SOLVERS)
Builds an EigenSolver using the linear solver package specified by solver_package.
Definition: eigen_solver.C:59
libMesh::DifferentiablePhysics::mass_residual
virtual bool mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on elem from elem_residual.
Definition: diff_physics.h:319
libMesh::ReferenceCounter::_enable_print_counter
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called.
Definition: reference_counter.h:141
libMesh::DifferentiablePhysics::side_time_derivative
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:172
libMesh::EigenTimeSolver::Matrix_A
The matrix associated with the spatial part of the operator.
Definition: eigen_time_solver.h:202
libMesh::TimeSolver::solution_history
std::unique_ptr< SolutionHistory > solution_history
A std::unique_ptr to a SolutionHistory object.
Definition: time_solver.h:264
libMesh::EigenTimeSolver::eigen_solver
std::unique_ptr< EigenSolver< Number > > eigen_solver
The EigenSolver object.
Definition: eigen_time_solver.h:155
libMesh::out
OStreamProxy out
libMesh::ImplicitSystem::get_matrix
const SparseMatrix< Number > & get_matrix(const std::string &mat_name) const
Definition: implicit_system.C:262
libMesh::TimeSolver::_linear_solver
std::unique_ptr< LinearSolver< Number > > _linear_solver
An implicit linear solver to use for adjoint problems.
Definition: time_solver.h:252
libMesh::Quality::name
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42