libMesh::EigenTimeSolver Class Reference

#include <eigen_time_solver.h>

Inheritance diagram for libMesh::EigenTimeSolver:

Public Types

typedef DifferentiableSystem sys_type
 
typedef TimeSolver Parent
 

Public Member Functions

 EigenTimeSolver (sys_type &s)
 
virtual ~EigenTimeSolver ()
 
virtual void init () libmesh_override
 
virtual void reinit () libmesh_override
 
virtual void solve () libmesh_override
 
virtual void advance_timestep () libmesh_override
 
Real error_order () const
 
virtual bool element_residual (bool get_jacobian, DiffContext &) libmesh_override
 
virtual bool side_residual (bool get_jacobian, DiffContext &) libmesh_override
 
virtual bool nonlocal_residual (bool get_jacobian, DiffContext &) libmesh_override
 
virtual Real du (const SystemNorm &) const libmesh_override
 
virtual bool is_steady () const libmesh_override
 
virtual void init_data ()
 
virtual void adjoint_advance_timestep ()
 
virtual void retrieve_timestep ()
 
virtual void before_timestep ()
 
const sys_typesystem () const
 
sys_typesystem ()
 
virtual UniquePtr< DiffSolver > & diff_solver ()
 
virtual UniquePtr< LinearSolver< Number > > & linear_solver ()
 
void set_solution_history (const SolutionHistory &_solution_history)
 
bool is_adjoint () const
 
void set_is_adjoint (bool _is_adjoint_value)
 

Static Public Member Functions

static std::string get_info ()
 
static void print_info (std::ostream &out=libMesh::out)
 
static unsigned int n_objects ()
 
static void enable_print_counter_info ()
 
static void disable_print_counter_info ()
 

Public Attributes

UniquePtr< EigenSolver< Number > > eigen_solver
 
Real tol
 
unsigned int maxits
 
unsigned int n_eigenpairs_to_compute
 
unsigned int n_basis_vectors_to_use
 
unsigned int n_converged_eigenpairs
 
unsigned int n_iterations_reqd
 
bool quiet
 
unsigned int reduce_deltat_on_diffsolver_failure
 

Protected Types

typedef bool(DifferentiablePhysics::* ResFuncType) (bool, DiffContext &)
 
typedef void(DiffContext::* ReinitFuncType) (Real)
 
typedef std::map< std::string, std::pair< unsigned int, unsigned int > > Counts
 

Protected Member Functions

void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

UniquePtr< DiffSolver_diff_solver
 
UniquePtr< LinearSolver< Number > > _linear_solver
 
sys_type_system
 
UniquePtr< SolutionHistorysolution_history
 

Static Protected Attributes

static Counts _counts
 
static Threads::atomic< unsigned int > _n_objects
 
static Threads::spin_mutex _mutex
 
static bool _enable_print_counter = true
 

Private Types

enum  NowAssembling { Matrix_A, Matrix_B, Invalid_Matrix }
 

Private Attributes

NowAssembling now_assembling
 

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

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 119 of file reference_counter.h.

The parent class

Definition at line 77 of file eigen_time_solver.h.

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

Definition at line 272 of file time_solver.h.

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.

The type of system

Definition at line 72 of file eigen_time_solver.h.

Member Enumeration Documentation

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.

Constructor & Destructor Documentation

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

Constructor. Requires a reference to the system to be solved.

Definition at line 32 of file eigen_time_solver.C.

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

33  : Parent(s),
35  tol(pow(TOLERANCE, 5./3.)),
36  maxits(1000),
41 {
42  libmesh_experimental();
43  eigen_solver->set_eigenproblem_type(GHEP);//or GNHEP
44  eigen_solver->set_position_of_spectrum(LARGEST_MAGNITUDE);
45 }
static const Real TOLERANCE
UniquePtr< EigenSolver< Number > > eigen_solver
double pow(double a, int b)
static UniquePtr< EigenSolver< T > > build(const Parallel::Communicator &comm_in LIBMESH_CAN_DEFAULT_TO_COMMWORLD, const SolverPackage solver_package=SLEPC_SOLVERS)
Definition: eigen_solver.C:37
unsigned int n_eigenpairs_to_compute
libMesh::EigenTimeSolver::~EigenTimeSolver ( )
virtual

Destructor.

Definition at line 47 of file eigen_time_solver.C.

48 {
49 }

Member Function Documentation

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 106 of file time_solver.C.

107 {
108 }
virtual void libMesh::EigenTimeSolver::advance_timestep ( )
inlinevirtual

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 {}
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 {}
virtual UniquePtr<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.

References libMesh::TimeSolver::_diff_solver.

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

181 { return _diff_solver; }
UniquePtr< DiffSolver > _diff_solver
Definition: time_solver.h:247
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited
virtual Real libMesh::EigenTimeSolver::du ( const SystemNorm ) const
inlinevirtual
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.; }
bool libMesh::EigenTimeSolver::element_residual ( bool  get_jacobian,
DiffContext context 
)
virtual

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

Implements libMesh::TimeSolver.

Definition at line 127 of file eigen_time_solver.C.

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.

Referenced by error_order().

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

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

Definition at line 101 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

102 {
103  _enable_print_counter = true;
104  return;
105 }
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.

References element_residual(), nonlocal_residual(), and side_residual().

118 { return 0.; }
std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

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

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 (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
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 185 of file reference_counter.h.

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

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

186 {
187  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
188  std::pair<unsigned int, unsigned int> & p = _counts[name];
189 
190  p.first++;
191 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
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 198 of file reference_counter.h.

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

Referenced by libMesh::ReferenceCounter::n_objects(), and libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

199 {
200  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
201  std::pair<unsigned int, unsigned int> & p = _counts[name];
202 
203  p.second++;
204 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::EigenTimeSolver::init ( )
virtual

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

Reimplemented from libMesh::TimeSolver.

Definition at line 56 of file eigen_time_solver.C.

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

57 {
58  // Add matrix "B" to _system if not already there.
59  // The user may have already added a matrix "B" before
60  // calling the System initialization. This would be
61  // necessary if e.g. the System originally started life
62  // with a different type of TimeSolver and only later
63  // had its TimeSolver changed to an EigenTimeSolver.
64  if (!_system.have_matrix("B"))
65  _system.add_matrix("B");
66 }
SparseMatrix< Number > & add_matrix(const std::string &mat_name)
sys_type & _system
Definition: time_solver.h:257
bool have_matrix(const std::string &mat_name) const
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 77 of file time_solver.C.

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().

78 {
79  this->diff_solver()->init();
80 
81  if (libMesh::on_command_line("--solver_system_names"))
82  this->linear_solver()->init((_system.name()+"_").c_str());
83  else
84  this->linear_solver()->init();
85 }
virtual UniquePtr< DiffSolver > & diff_solver()
Definition: time_solver.h:181
const std::string & name() const
Definition: system.h:1996
sys_type & _system
Definition: time_solver.h:257
bool on_command_line(const std::string &arg)
Definition: libmesh.C:921
virtual UniquePtr< LinearSolver< Number > > & linear_solver()
Definition: time_solver.h:186
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.

References libMesh::TimeSolver::_is_adjoint.

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

233  { return _is_adjoint; }
virtual bool libMesh::EigenTimeSolver::is_steady ( ) const
inlinevirtual

This is effectively a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 149 of file eigen_time_solver.h.

149 { return true; }
virtual UniquePtr<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.

References libMesh::TimeSolver::_linear_solver.

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

186 { return _linear_solver; }
UniquePtr< LinearSolver< Number > > _linear_solver
Definition: time_solver.h:252
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited
bool libMesh::EigenTimeSolver::nonlocal_residual ( bool  get_jacobian,
DiffContext context 
)
virtual

Forms the jacobian of the nonlocal terms.

Implements libMesh::TimeSolver.

Definition at line 221 of file eigen_time_solver.C.

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.

Referenced by error_order().

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

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

Definition at line 88 of file reference_counter.C.

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

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

89 {
91  out_stream << ReferenceCounter::get_info();
92 }
static std::string get_info()
void libMesh::EigenTimeSolver::reinit ( )
virtual

The reinitialization function. This method is used after changes in the mesh

Reimplemented from libMesh::TimeSolver.

Definition at line 51 of file eigen_time_solver.C.

52 {
53  // empty...
54 }
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 110 of file time_solver.C.

111 {
112 }
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.

References libMesh::TimeSolver::_is_adjoint.

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

240  { _is_adjoint = _is_adjoint_value; }
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 97 of file time_solver.C.

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

98 {
99  solution_history = _solution_history.clone();
100 }
UniquePtr< SolutionHistory > solution_history
Definition: time_solver.h:264
bool libMesh::EigenTimeSolver::side_residual ( bool  get_jacobian,
DiffContext context 
)
virtual

Forms the jacobian of the boundary terms.

Implements libMesh::TimeSolver.

Definition at line 174 of file eigen_time_solver.C.

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().

Referenced by error_order().

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

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

Reimplemented from libMesh::TimeSolver.

Definition at line 68 of file eigen_time_solver.C.

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.

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

References libMesh::TimeSolver::_system.

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

171 { return _system; }
sys_type & _system
Definition: time_solver.h:257
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.

References libMesh::TimeSolver::_system.

176 { return _system; }
sys_type & _system
Definition: time_solver.h:257

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
UniquePtr<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().

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 143 of file reference_counter.h.

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

UniquePtr<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().

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 137 of file reference_counter.h.

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 132 of file reference_counter.h.

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

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().

UniquePtr<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().

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().

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().

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().

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().

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().

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().

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().

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().

UniquePtr<SolutionHistory> libMesh::TimeSolver::solution_history
protectedinherited

A UniquePtr 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().

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: