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

This class defines a Newmark time integrator for second order (in time) DifferentiableSystems. More...

#include <newmark_solver.h>

Inheritance diagram for libMesh::NewmarkSolver:
[legend]

Public Types

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

Public Member Functions

 NewmarkSolver (sys_type &s)
 Constructor. More...
 
virtual ~NewmarkSolver ()
 Destructor. More...
 
virtual void advance_timestep () override
 This method advances the solution to the next timestep, after a solve() has been performed. More...
 
virtual void adjoint_advance_timestep () override
 This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. More...
 
virtual void compute_initial_accel ()
 This method uses the specified initial displacement and velocity to compute the initial acceleration \(a_0\). More...
 
void project_initial_accel (FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr)
 Specify non-zero initial acceleration. More...
 
void set_initial_accel_avail (bool initial_accel_set)
 Allow the user to (re)set whether the initial acceleration is available. More...
 
virtual Real error_order () const override
 Error convergence order: 2 for \(\gamma=0.5\), 1 otherwise. More...
 
virtual void solve () override
 This method solves for the solution at the next timestep. More...
 
virtual bool element_residual (bool request_jacobian, DiffContext &) override
 This method uses the DifferentiablePhysics' element_time_derivative() and element_constraint() to build a full residual on an element. More...
 
virtual bool side_residual (bool request_jacobian, DiffContext &) override
 This method uses the DifferentiablePhysics' side_time_derivative() and side_constraint() to build a full residual on an element's side. More...
 
virtual bool nonlocal_residual (bool request_jacobian, DiffContext &) override
 This method uses the DifferentiablePhysics' nonlocal_time_derivative() and nonlocal_constraint() to build a full residual for non-local terms. More...
 
void set_beta (Real beta)
 Setter for \( \beta \). More...
 
void set_gamma (Real gamma)
 Setter for \( \gamma \). More...
 
virtual unsigned int time_order () const override
 
virtual void init () override
 The initialization function. More...
 
virtual void init_data () override
 The data initialization function. More...
 
virtual void reinit () override
 The reinitialization function. More...
 
virtual void retrieve_timestep () override
 This method retrieves all the stored solutions at the current system.time. More...
 
void project_initial_rate (FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr)
 Specify non-zero initial velocity. More...
 
Number old_solution_rate (const dof_id_type global_dof_number) const
 
Number old_solution_accel (const dof_id_type global_dof_number) const
 
virtual void init_adjoints () override
 Add adjoint vectors and old_adjoint_vectors as per the indices of QoISet. More...
 
void update ()
 
virtual std::pair< unsigned int, Realadjoint_solve (const QoISet &qoi_indices) override
 This method solves for the adjoint solution at the next adjoint timestep (or a steady state adjoint solve) More...
 
virtual void integrate_qoi_timestep () override
 A method to integrate the system::QoI functionals. More...
 
virtual void integrate_adjoint_sensitivity (const QoISet &qois, const ParameterVector &parameter_vector, SensitivityData &sensitivities) override
 A method to integrate the adjoint sensitivity w.r.t a given parameter vector. More...
 
virtual void integrate_adjoint_refinement_error_estimate (AdjointRefinementEstimator &, ErrorVector &) override
 A method to compute the adjoint refinement error estimate at the current timestep. More...
 
Number old_nonlinear_solution (const dof_id_type global_dof_number) const
 
virtual Real du (const SystemNorm &norm) const override
 Computes the size of ||u^{n+1} - u^{n}|| in some norm. More...
 
virtual bool is_steady () const override
 This is not a steady-state solver. More...
 
void set_first_adjoint_step (bool first_adjoint_step_setting)
 A setter for the first_adjoint_step boolean. More...
 
void set_first_solve (bool first_solve_setting)
 
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...
 
SolutionHistoryget_solution_history ()
 A getter function that returns a reference to the solution history object owned by TimeSolver. 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...
 
virtual Real last_completed_timestep_size ()
 Returns system.deltat if fixed timestep solver is used, the complete timestep size (sum of all substeps) if the adaptive time solver is used. 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_stream=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::shared_ptr< NumericVector< Number > > old_local_nonlinear_solution
 Serial vector of _system.get_vector("_old_nonlinear_solution") This is a shared_ptr so that it can be shared between different derived class instances, as in e.g. 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

virtual bool _general_residual (bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit)
 This method is the underlying implementation of the public residual methods. More...
 
void increment_constructor_count (const std::string &name) noexcept
 Increments the construction counter. More...
 
void increment_destructor_count (const std::string &name) noexcept
 Increments the destruction counter. More...
 

Protected Attributes

Real _beta
 The value for the \(\beta\) to employ. More...
 
Real _gamma
 The value for \(\gamma\) to employ. More...
 
bool _is_accel_solve
 Need to be able to indicate to _general_residual if we are doing an acceleration solve or not. More...
 
bool _initial_accel_set
 This method requires an initial acceleration. More...
 
std::unique_ptr< NumericVector< Number > > _old_local_solution_rate
 Serial vector of previous time step velocity \( \dot{u}_n \). More...
 
std::unique_ptr< NumericVector< Number > > _old_local_solution_accel
 Serial vector of previous time step acceleration \( \ddot{u}_n \). More...
 
bool first_solve
 A bool that will be true the first time solve() is called, and false thereafter. More...
 
bool first_adjoint_step
 A bool that will be true the first time adjoint_advance_timestep() is called, (when the primal solution is to be used to set adjoint boundary conditions) and false thereafter. More...
 
std::vector< std::unique_ptr< NumericVector< Number > > > old_adjoints
 A vector of pointers to vectors holding the adjoint solution at the last time step. More...
 
Real last_step_deltat
 We will need to move the system.time around to ensure that residuals are built with the right deltat and the right time. More...
 
Real next_step_deltat
 
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...
 
Real last_deltat
 The deltat for the last completed timestep before the current one. 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...
 

Detailed Description

This class defines a Newmark time integrator for second order (in time) DifferentiableSystems.

There are two parameters \(\gamma\) and \(\beta\) (defaulting to 0.5 and 0.25, respectively).

Note
Projectors are included for the initial velocity and acceleration; the initial displacement can be set by calling System::project_solution().

This class is part of the new DifferentiableSystem framework, which is still experimental. Users of this framework should beware of bugs and future API changes.

Author
Paul T. Bauman
Date
2015

Definition at line 46 of file newmark_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 119 of file reference_counter.h.

◆ Parent

The parent class.

Definition at line 52 of file newmark_solver.h.

◆ ReinitFuncType

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

Definition at line 327 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 325 of file time_solver.h.

◆ sys_type

The type of system.

Definition at line 69 of file time_solver.h.

Constructor & Destructor Documentation

◆ NewmarkSolver()

libMesh::NewmarkSolver::NewmarkSolver ( sys_type s)
explicit

Constructor.

Requires a reference to the system to be solved.

Definition at line 27 of file newmark_solver.C.

29  _beta(0.25),
30  _gamma(0.5),
31  _is_accel_solve(false),
32  _initial_accel_set(false)
33 {}
Real _gamma
The value for to employ.
bool _initial_accel_set
This method requires an initial acceleration.
SecondOrderUnsteadySolver(sys_type &s)
Constructor.
bool _is_accel_solve
Need to be able to indicate to _general_residual if we are doing an acceleration solve or not...
Real _beta
The value for the to employ.

◆ ~NewmarkSolver()

libMesh::NewmarkSolver::~NewmarkSolver ( )
virtualdefault

Destructor.

Member Function Documentation

◆ _general_residual()

bool libMesh::NewmarkSolver::_general_residual ( bool  request_jacobian,
DiffContext context,
ResFuncType  mass,
ResFuncType  damping,
ResFuncType  time_deriv,
ResFuncType  constraint,
ReinitFuncType  reinit 
)
protectedvirtual

This method is the underlying implementation of the public residual methods.

Definition at line 203 of file newmark_solver.C.

References _beta, _gamma, _is_accel_solve, libMesh::TimeSolver::_system, libMesh::DenseVector< T >::add(), libMesh::DifferentiableSystem::deltat, libMesh::DiffContext::elem_solution_accel_derivative, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_fixed_solution(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_solution(), libMesh::DiffContext::get_elem_solution_accel(), libMesh::DiffContext::get_elem_solution_rate(), libMesh::DifferentiableSystem::get_physics(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::SecondOrderUnsteadySolver::old_solution_accel(), libMesh::SecondOrderUnsteadySolver::old_solution_rate(), libMesh::Real, libMesh::DenseVector< T >::size(), libMesh::DenseMatrix< T >::swap(), and libMesh::System::use_fixed_solution.

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

210 {
211  unsigned int n_dofs = context.get_elem_solution().size();
212 
213  // We might need to save the old jacobian in case one of our physics
214  // terms later is unable to update it analytically.
215  DenseMatrix<Number> old_elem_jacobian(n_dofs, n_dofs);
216 
217  // Local velocity at old time step
218  DenseVector<Number> old_elem_solution_rate(n_dofs);
219  for (unsigned int i=0; i != n_dofs; ++i)
220  old_elem_solution_rate(i) =
221  old_solution_rate(context.get_dof_indices()[i]);
222 
223  // The user is computing the initial acceleration
224  // So upstream we've swapped _system.solution and _old_local_solution_accel
225  // So we need to give the context the correct entries since we're solving for
226  // acceleration here.
227  if (_is_accel_solve)
228  {
229  // System._solution is actually the acceleration right now so we need
230  // to reset the elem_solution to the right thing, which in this case
231  // is the initial guess for displacement, which got swapped into
232  // _old_solution_accel vector
233  DenseVector<Number> old_elem_solution(n_dofs);
234  for (unsigned int i=0; i != n_dofs; ++i)
235  old_elem_solution(i) =
236  old_solution_accel(context.get_dof_indices()[i]);
237 
238  context.elem_solution_derivative = 0.0;
239  context.elem_solution_rate_derivative = 0.0;
240  context.elem_solution_accel_derivative = 1.0;
241 
242  // Acceleration is currently the unknown so it's already sitting
243  // in elem_solution() thanks to FEMContext::pre_fe_reinit
244  context.get_elem_solution_accel() = context.get_elem_solution();
245 
246  // Now reset elem_solution() to what the user is expecting
247  context.get_elem_solution() = old_elem_solution;
248 
249  context.get_elem_solution_rate() = old_elem_solution_rate;
250 
251  // The user's Jacobians will be targeting derivatives w.r.t. u_{n+1}.
252  // Although the vast majority of cases will have the correct analytic
253  // Jacobians in this iteration, since we reset elem_solution_derivative*,
254  // if there are coupled/overlapping problems, there could be
255  // mismatches in the Jacobian. So we force finite differencing for
256  // the first iteration.
257  request_jacobian = false;
258  }
259  // Otherwise, the unknowns are the displacements and everything is straight
260  // forward and is what you think it is
261  else
262  {
263  if (request_jacobian)
264  old_elem_jacobian.swap(context.get_elem_jacobian());
265 
266  // Local displacement at old timestep
267  DenseVector<Number> old_elem_solution(n_dofs);
268  for (unsigned int i=0; i != n_dofs; ++i)
269  old_elem_solution(i) =
270  old_nonlinear_solution(context.get_dof_indices()[i]);
271 
272  // Local acceleration at old time step
273  DenseVector<Number> old_elem_solution_accel(n_dofs);
274  for (unsigned int i=0; i != n_dofs; ++i)
275  old_elem_solution_accel(i) =
276  old_solution_accel(context.get_dof_indices()[i]);
277 
278  // Convenience
280 
281  context.elem_solution_derivative = 1.0;
282 
283  // Local velocity at current time step
284  // v_{n+1} = gamma/(beta*Delta t)*(x_{n+1}-x_n)
285  // + (1-(gamma/beta))*v_n
286  // + (1-gamma/(2*beta))*(Delta t)*a_n
287  context.elem_solution_rate_derivative = (_gamma/(_beta*dt));
288 
289  context.get_elem_solution_rate() = context.get_elem_solution();
290  context.get_elem_solution_rate() -= old_elem_solution;
291  context.get_elem_solution_rate() *= context.elem_solution_rate_derivative;
292  context.get_elem_solution_rate().add( (1.0-_gamma/_beta), old_elem_solution_rate);
293  context.get_elem_solution_rate().add( (1.0-_gamma/(2.0*_beta))*dt, old_elem_solution_accel);
294 
295 
296 
297  // Local acceleration at current time step
298  // a_{n+1} = (1/(beta*(Delta t)^2))*(x_{n+1}-x_n)
299  // - 1/(beta*Delta t)*v_n
300  // - (1/(2*beta)-1)*a_n
301  context.elem_solution_accel_derivative = 1.0/(_beta*dt*dt);
302 
303  context.get_elem_solution_accel() = context.get_elem_solution();
304  context.get_elem_solution_accel() -= old_elem_solution;
305  context.get_elem_solution_accel() *= context.elem_solution_accel_derivative;
306  context.get_elem_solution_accel().add(-1.0/(_beta*dt), old_elem_solution_rate);
307  context.get_elem_solution_accel().add(-(1.0/(2.0*_beta)-1.0), old_elem_solution_accel);
308 
309  // Move the mesh into place first if necessary, set t = t_{n+1}
310  (context.*reinit_func)(1.);
311  }
312 
313  // If a fixed solution is requested, we'll use x_{n+1}
315  context.get_elem_fixed_solution() = context.get_elem_solution();
316 
317  // Get the time derivative at t_{n+1}, F(u_{n+1})
318  bool jacobian_computed = (_system.get_physics()->*time_deriv)(request_jacobian, context);
319 
320  // Damping at t_{n+1}, C(u_{n+1})
321  jacobian_computed = (_system.get_physics()->*damping)(jacobian_computed, context) &&
322  jacobian_computed;
323 
324  // Mass at t_{n+1}, M(u_{n+1})
325  jacobian_computed = (_system.get_physics()->*mass)(jacobian_computed, context) &&
326  jacobian_computed;
327 
328  // Add the constraint term
329  jacobian_computed = (_system.get_physics()->*constraint)(jacobian_computed, context) &&
330  jacobian_computed;
331 
332  // Add back (or restore) the old jacobian
333  if (request_jacobian)
334  {
335  if (jacobian_computed)
336  context.get_elem_jacobian() += old_elem_jacobian;
337  else
338  context.get_elem_jacobian().swap(old_elem_jacobian);
339  }
340 
341  return jacobian_computed;
342 }
Real _gamma
The value for to employ.
Number old_solution_accel(const dof_id_type global_dof_number) const
Number old_nonlinear_solution(const dof_id_type global_dof_number) const
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
Definition: diff_system.h:278
bool use_fixed_solution
A boolean to be set to true by systems using elem_fixed_solution, for optional use by e...
Definition: system.h:1543
Number old_solution_rate(const dof_id_type global_dof_number) const
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:181
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
bool _is_accel_solve
Need to be able to indicate to _general_residual if we are doing an acceleration solve or not...
Real _beta
The value for the to employ.

◆ adjoint_advance_timestep()

void libMesh::NewmarkSolver::adjoint_advance_timestep ( )
overridevirtual

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 from libMesh::UnsteadySolver.

Definition at line 102 of file newmark_solver.C.

103 {
104  libmesh_not_implemented();
105 }

◆ adjoint_solve()

std::pair< unsigned int, Real > libMesh::UnsteadySolver::adjoint_solve ( const QoISet qoi_indices)
overridevirtualinherited

This method solves for the adjoint solution at the next adjoint timestep (or a steady state adjoint solve)

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.

Definition at line 228 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, and libMesh::TimeSolver::last_deltat.

229 {
230  std::pair<unsigned int, Real> adjoint_output = _system.ImplicitSystem::adjoint_solve(qoi_indices);
231 
232  // Record the deltat we used for this adjoint timestep. This was determined completely
233  // by SolutionHistory::retrieve methods. The adjoint_solve methods should never change deltat.
235 
236  return adjoint_output;
237 }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
Definition: diff_system.h:278
Real last_deltat
The deltat for the last completed timestep before the current one.
Definition: time_solver.h:332

◆ advance_timestep()

void libMesh::NewmarkSolver::advance_timestep ( )
overridevirtual

This method advances the solution to the next timestep, after a solve() has been performed.

Often this will be done after every UnsteadySolver::solve(), but adaptive mesh refinement and/or adaptive time step selection may require some solve() steps to be repeated.

Reimplemented from libMesh::UnsteadySolver.

Definition at line 44 of file newmark_solver.C.

References _beta, _gamma, libMesh::SecondOrderUnsteadySolver::_old_local_solution_accel, libMesh::SecondOrderUnsteadySolver::_old_local_solution_rate, libMesh::TimeSolver::_system, libMesh::NumericVector< T >::add(), libMesh::UnsteadySolver::advance_timestep(), libMesh::NumericVector< T >::clone(), libMesh::DifferentiableSystem::deltat, libMesh::UnsteadySolver::first_solve, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::SecondOrderUnsteadySolver::old_solution_accel(), libMesh::SecondOrderUnsteadySolver::old_solution_rate(), and libMesh::System::solution.

45 {
46  // We need to update velocity and acceleration before
47  // we update the nonlinear solution (displacement) and
48  // delta_t
49 
51  _system.get_vector("_old_solution_rate");
52 
54  _system.get_vector("_old_solution_accel");
55 
56  if (!first_solve)
57  {
58  NumericVector<Number> & old_nonlinear_soln =
59  _system.get_vector("_old_nonlinear_solution");
60 
61  NumericVector<Number> & nonlinear_solution =
62  *(_system.solution);
63 
64  // We need to cache the new solution_rate before updating the old_solution_rate
65  // so we can update acceleration with the proper old_solution_rate
66  // v_{n+1} = gamma/(beta*Delta t)*(x_{n+1}-x_n)
67  // - ((gamma/beta)-1)*v_n
68  // - (gamma/(2*beta)-1)*(Delta t)*a_n
69  std::unique_ptr<NumericVector<Number>> new_solution_rate = nonlinear_solution.clone();
70  (*new_solution_rate) -= old_nonlinear_soln;
71  (*new_solution_rate) *= (_gamma/(_beta*_system.deltat));
72  new_solution_rate->add( (1.0-_gamma/_beta), old_solution_rate );
73  new_solution_rate->add( (1.0-_gamma/(2.0*_beta))*_system.deltat, old_solution_accel );
74 
75  // a_{n+1} = (1/(beta*(Delta t)^2))*(x_{n+1}-x_n)
76  // - 1/(beta*Delta t)*v_n
77  // - (1-1/(2*beta))*a_n
78  std::unique_ptr<NumericVector<Number>> new_solution_accel = old_solution_accel.clone();
79  (*new_solution_accel) *= -(1.0/(2.0*_beta)-1.0);
80  new_solution_accel->add( -1.0/(_beta*_system.deltat), old_solution_rate );
81  new_solution_accel->add( 1.0/(_beta*_system.deltat*_system.deltat), nonlinear_solution );
82  new_solution_accel->add( -1.0/(_beta*_system.deltat*_system.deltat), old_nonlinear_soln );
83 
84  // Now update old_solution_rate
85  old_solution_rate = (*new_solution_rate);
86  old_solution_accel = (*new_solution_accel);
87  }
88 
89  // Localize updated vectors
90  old_solution_rate.localize
93 
94  old_solution_accel.localize
97 
98  // Now we can finish advancing the timestep
100 }
Real _gamma
The value for to employ.
Number old_solution_accel(const dof_id_type global_dof_number) const
std::unique_ptr< NumericVector< Number > > _old_local_solution_accel
Serial vector of previous time step acceleration .
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
Definition: diff_system.h:278
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1573
Number old_solution_rate(const dof_id_type global_dof_number) const
virtual void advance_timestep() override
This method advances the solution to the next timestep, after a solve() has been performed.
std::unique_ptr< NumericVector< Number > > _old_local_solution_rate
Serial vector of previous time step velocity .
const DofMap & get_dof_map() const
Definition: system.h:2293
template class LIBMESH_EXPORT NumericVector< Number >
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:511
Real _beta
The value for the to employ.
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918
bool first_solve
A bool that will be true the first time solve() is called, and false thereafter.

◆ 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 205 of file time_solver.h.

205 {}

◆ compute_initial_accel()

void libMesh::NewmarkSolver::compute_initial_accel ( )
virtual

This method uses the specified initial displacement and velocity to compute the initial acceleration \(a_0\).

Definition at line 107 of file newmark_solver.C.

References libMesh::TimeSolver::_diff_solver, _is_accel_solve, libMesh::TimeSolver::_system, libMesh::System::get_vector(), set_initial_accel_avail(), libMesh::System::solution, and libMesh::System::update().

Referenced by main().

108 {
109  // We need to compute the initial acceleration based off of
110  // the initial position and velocity and, thus, acceleration
111  // is the unknown in diff_solver and not the displacement. So,
112  // We swap solution and acceleration. NewmarkSolver::_general_residual
113  // will check _is_accel_solve and provide the correct
114  // values to the FEMContext assuming this swap was made.
115  this->_is_accel_solve = true;
116 
117  //solution->accel, accel->solution
118  _system.solution->swap(_system.get_vector("_old_solution_accel"));
119  _system.update();
120 
121  this->_diff_solver->solve();
122 
123  // solution->solution, accel->accel
124  _system.solution->swap(_system.get_vector("_old_solution_accel"));
125  _system.update();
126 
127  // We're done, so no longer doing an acceleration solve
128  this->_is_accel_solve = false;
129 
130  this->set_initial_accel_avail(true);
131 }
std::unique_ptr< DiffSolver > _diff_solver
An implicit linear or nonlinear solver to use at each timestep.
Definition: time_solver.h:302
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1573
void set_initial_accel_avail(bool initial_accel_set)
Allow the user to (re)set whether the initial acceleration is available.
virtual void update()
Update the local values to reflect the solution on neighboring processors.
Definition: system.C:493
bool _is_accel_solve
Need to be able to indicate to _general_residual if we are doing an acceleration solve or not...
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918

◆ 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 220 of file time_solver.h.

References libMesh::TimeSolver::_diff_solver.

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

220 { return _diff_solver; }
std::unique_ptr< DiffSolver > _diff_solver
An implicit linear or nonlinear solver to use at each timestep.
Definition: time_solver.h:302

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

101 {
102  _enable_print_counter = false;
103  return;
104 }
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ du()

Real libMesh::UnsteadySolver::du ( const SystemNorm norm) const
overridevirtualinherited

Computes the size of ||u^{n+1} - u^{n}|| in some norm.

Note
While you can always call this function, its result may or may not be very meaningful. For example, if you call this function right after calling advance_timestep() then you'll get a result of zero since old_nonlinear_solution is set equal to nonlinear_solution in this function.

Implements libMesh::TimeSolver.

Definition at line 348 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::calculate_norm(), libMesh::System::get_vector(), libMesh::TensorTools::norm(), and libMesh::System::solution.

349 {
350 
351  std::unique_ptr<NumericVector<Number>> solution_copy =
352  _system.solution->clone();
353 
354  solution_copy->add(-1., _system.get_vector("_old_nonlinear_solution"));
355 
356  solution_copy->close();
357 
358  return _system.calculate_norm(*solution_copy, norm);
359 }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
std::unique_ptr< NumericVector< Number > > solution
Data structure to hold solution values.
Definition: system.h:1573
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
Definition: system.C:1672
auto norm(const T &a) -> decltype(std::abs(a))
Definition: tensor_tools.h:74
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918

◆ element_residual()

bool libMesh::NewmarkSolver::element_residual ( bool  request_jacobian,
DiffContext context 
)
overridevirtual

This method uses the DifferentiablePhysics' element_time_derivative() and element_constraint() to build a full residual on an element.

What combination it uses will depend on theta.

Implements libMesh::TimeSolver.

Definition at line 161 of file newmark_solver.C.

References libMesh::DifferentiablePhysics::_eulerian_time_deriv(), _general_residual(), libMesh::DifferentiablePhysics::damping_residual(), libMesh::DiffContext::elem_reinit(), libMesh::DifferentiablePhysics::element_constraint(), and libMesh::DifferentiablePhysics::mass_residual().

163 {
164  return this->_general_residual(request_jacobian,
165  context,
171 }
virtual bool damping_residual(bool request_jacobian, DiffContext &)
Subtracts a damping vector contribution on elem from elem_residual.
Definition: diff_physics.h:360
virtual void elem_reinit(Real)
Gives derived classes the opportunity to reinitialize data (FE objects in FEMSystem, for example) needed for an interior integration at a new point within a timestep.
Definition: diff_context.h:77
bool _eulerian_time_deriv(bool request_jacobian, DiffContext &)
This method simply combines element_time_derivative() and eulerian_residual(), which makes its addres...
Definition: diff_physics.C:96
virtual bool mass_residual(bool request_jacobian, DiffContext &)
Subtracts a mass vector contribution on elem from elem_residual.
Definition: diff_physics.h:302
virtual bool _general_residual(bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit)
This method is the underlying implementation of the public residual methods.
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on elem to elem_residual.
Definition: diff_physics.h:144

◆ 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 94 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

95 {
96  _enable_print_counter = true;
97  return;
98 }
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ error_order()

Real libMesh::NewmarkSolver::error_order ( ) const
overridevirtual

Error convergence order: 2 for \(\gamma=0.5\), 1 otherwise.

Implements libMesh::UnsteadySolver.

Definition at line 37 of file newmark_solver.C.

References _gamma.

38 {
39  if (_gamma == 0.5)
40  return 2.;
41  return 1.;
42 }
Real _gamma
The value for to employ.

◆ 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.

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 (const auto & [name, cd] : _counts)
59  oss << "| " << name << " reference count information:\n"
60  << "| Creations: " << cd.first << '\n'
61  << "| Destructions: " << cd.second << '\n';
62 
63  oss << " ---------------------------------------------------------------------------- \n";
64 
65  return oss.str();
66 
67 #else
68 
69  return "";
70 
71 #endif
72 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
static Counts _counts
Actually holds the data.

◆ get_solution_history()

SolutionHistory & libMesh::TimeSolver::get_solution_history ( )
inherited

A getter function that returns a reference to the solution history object owned by TimeSolver.

Definition at line 124 of file time_solver.C.

References libMesh::TimeSolver::solution_history.

Referenced by libMesh::AdaptiveTimeSolver::init().

125 {
126  return *solution_history;
127 }
std::unique_ptr< SolutionHistory > solution_history
A std::unique_ptr to a SolutionHistory object.
Definition: time_solver.h:319

◆ increment_constructor_count()

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

Increments the construction counter.

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

Definition at line 183 of file reference_counter.h.

References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

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

184 {
185  libmesh_try
186  {
187  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
188  std::pair<unsigned int, unsigned int> & p = _counts[name];
189  p.first++;
190  }
191  libmesh_catch (...)
192  {
193  auto stream = libMesh::err.get();
194  stream->exceptions(stream->goodbit); // stream must not throw
195  libMesh::err << "Encountered unrecoverable error while calling "
196  << "ReferenceCounter::increment_constructor_count() "
197  << "for a(n) " << name << " object." << std::endl;
198  std::terminate();
199  }
200 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
OStreamProxy err
static Counts _counts
Actually holds the data.
streamT * get()
Rather than implement every ostream/ios/ios_base function, we&#39;ll be lazy and make esoteric uses go th...
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30

◆ increment_destructor_count()

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

Increments the destruction counter.

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

Definition at line 207 of file reference_counter.h.

References libMesh::err, libMesh::BasicOStreamProxy< charT, traits >::get(), libMesh::Quality::name(), and libMesh::Threads::spin_mtx.

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

208 {
209  libmesh_try
210  {
211  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
212  std::pair<unsigned int, unsigned int> & p = _counts[name];
213  p.second++;
214  }
215  libmesh_catch (...)
216  {
217  auto stream = libMesh::err.get();
218  stream->exceptions(stream->goodbit); // stream must not throw
219  libMesh::err << "Encountered unrecoverable error while calling "
220  << "ReferenceCounter::increment_destructor_count() "
221  << "for a(n) " << name << " object." << std::endl;
222  std::terminate();
223  }
224 }
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:42
OStreamProxy err
static Counts _counts
Actually holds the data.
streamT * get()
Rather than implement every ostream/ios/ios_base function, we&#39;ll be lazy and make esoteric uses go th...
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
Definition: threads.C:30

◆ init()

void libMesh::SecondOrderUnsteadySolver::init ( )
overridevirtualinherited

The initialization function.

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

Reimplemented from libMesh::UnsteadySolver.

Definition at line 34 of file second_order_unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::add_vector(), and libMesh::UnsteadySolver::init().

35 {
37 
38  _system.add_vector("_old_solution_rate");
39  _system.add_vector("_old_solution_accel");
40 }
NumericVector< Number > & add_vector(std::string_view vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
Definition: system.C:751
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
virtual void init() override
The initialization function.

◆ init_adjoints()

void libMesh::UnsteadySolver::init_adjoints ( )
overridevirtualinherited

Add adjoint vectors and old_adjoint_vectors as per the indices of QoISet.

Reimplemented from libMesh::TimeSolver.

Definition at line 68 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::add_vector(), libMesh::GHOSTED, libMesh::TimeSolver::init_adjoints(), libMesh::make_range(), and libMesh::System::n_qois().

69 {
71 
72  // Add old adjoint solutions
73  // To keep the number of vectors consistent between the primal and adjoint
74  // time loops, we will also add the adjoint rhs vector during initialization
75  for(auto i : make_range(_system.n_qois()))
76  {
77  std::string old_adjoint_solution_name = "_old_adjoint_solution";
78  old_adjoint_solution_name+= std::to_string(i);
79  _system.add_vector(old_adjoint_solution_name, false, GHOSTED);
80 
81  std::string adjoint_rhs_name = "adjoint_rhs";
82  adjoint_rhs_name+= std::to_string(i);
83  _system.add_vector(adjoint_rhs_name, false, GHOSTED);
84  }
85 
86 }
unsigned int n_qois() const
Number of currently active quantities of interest.
Definition: system.h:2516
NumericVector< Number > & add_vector(std::string_view vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Adds the additional vector vec_name to this system.
Definition: system.C:751
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134
virtual void init_adjoints()
Initialize any adjoint related data structures, based on the number of qois.
Definition: time_solver.C:83

◆ init_data()

void libMesh::SecondOrderUnsteadySolver::init_data ( )
overridevirtualinherited

The data initialization function.

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

Reimplemented from libMesh::UnsteadySolver.

Definition at line 42 of file second_order_unsteady_solver.C.

References libMesh::SecondOrderUnsteadySolver::_old_local_solution_accel, libMesh::SecondOrderUnsteadySolver::_old_local_solution_rate, libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::GHOSTED, libMesh::UnsteadySolver::init_data(), libMesh::System::n_dofs(), libMesh::System::n_local_dofs(), and libMesh::SERIAL.

43 {
45 
46 #ifdef LIBMESH_ENABLE_GHOSTED
49  GHOSTED);
50 
53  GHOSTED);
54 #else
57 #endif
58 }
virtual void init_data() override
The data initialization function.
std::unique_ptr< NumericVector< Number > > _old_local_solution_accel
Serial vector of previous time step acceleration .
dof_id_type n_local_dofs() const
Definition: system.C:150
dof_id_type n_dofs() const
Definition: system.C:113
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
std::unique_ptr< NumericVector< Number > > _old_local_solution_rate
Serial vector of previous time step velocity .
const DofMap & get_dof_map() const
Definition: system.h:2293
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:511

◆ integrate_adjoint_refinement_error_estimate()

void libMesh::UnsteadySolver::integrate_adjoint_refinement_error_estimate ( AdjointRefinementEstimator ,
ErrorVector  
)
overridevirtualinherited

A method to compute the adjoint refinement error estimate at the current timestep.

int_{tstep_start}^{tstep_end} R(u^h,z) dt The user provides an initialized ARefEE object. Fills in an ErrorVector that contains the weighted sum of errors from all the QoIs and can be used to guide AMR. CURRENTLY ONLY SUPPORTED for Backward Euler.

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::EulerSolver, libMesh::FirstOrderUnsteadySolver, libMesh::TwostepTimeSolver, libMesh::AdaptiveTimeSolver, and libMesh::Euler2Solver.

Definition at line 331 of file unsteady_solver.C.

332 {
333  libmesh_not_implemented();
334 }

◆ integrate_adjoint_sensitivity()

void libMesh::UnsteadySolver::integrate_adjoint_sensitivity ( const QoISet qois,
const ParameterVector parameter_vector,
SensitivityData sensitivities 
)
overridevirtualinherited

A method to integrate the adjoint sensitivity w.r.t a given parameter vector.

int_{tstep_start}^{tstep_end} dQ/dp dt = int_{tstep_start}^{tstep_end} ( / p) - ( R (u,z) / p ) dt The trapezoidal rule is used to numerically integrate the timestep.

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.

Definition at line 280 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::ImplicitSystem::adjoint_qoi_parameter_sensitivity(), libMesh::DifferentiableSystem::deltat, libMesh::make_range(), libMesh::Real, libMesh::System::remove_vector(), libMesh::UnsteadySolver::retrieve_timestep(), libMesh::ParameterVector::size(), libMesh::QoISet::size(), and libMesh::System::time.

281 {
282  // CURRENTLY using the trapezoidal rule to integrate each timestep
283  // (f(t_j) + f(t_j+1))/2 (t_j+1 - t_j)
284  // Fix me: This function needs to be moved to the EulerSolver classes like the
285  // other integrate_timestep functions, and use an integration rule consistent with
286  // the theta method used for the time integration.
287 
288  // Get t_j
289  Real time_left = _system.time;
290 
291  // Left side sensitivities to hold f(t_j)
292  SensitivityData sensitivities_left(qois, _system, parameter_vector);
293 
294  // Get f(t_j)
295  _system.adjoint_qoi_parameter_sensitivity(qois, parameter_vector, sensitivities_left);
296 
297  // Advance to t_j+1
299 
300  // Get t_j+1
301  Real time_right = _system.time;
302 
303  // Right side sensitivities f(t_j+1)
304  SensitivityData sensitivities_right(qois, _system, parameter_vector);
305 
306  // Remove the sensitivity rhs vector from system since we did not write it to file and it cannot be retrieved
307  _system.remove_vector("sensitivity_rhs0");
308 
309  // Retrieve the primal and adjoint solutions at the current timestep
311 
312  // Get f(t_j+1)
313  _system.adjoint_qoi_parameter_sensitivity(qois, parameter_vector, sensitivities_right);
314 
315  // Remove the sensitivity rhs vector from system since we did not write it to file and it cannot be retrieved
316  _system.remove_vector("sensitivity_rhs0");
317 
318  // Get the contributions for each sensitivity from this timestep
319  const auto pv_size = parameter_vector.size();
320  for (auto i : make_range(qois.size(_system)))
321  for (auto j : make_range(pv_size))
322  sensitivities[i][j] = ( (sensitivities_left[i][j] + sensitivities_right[i][j])/2. )*(time_right - time_left);
323 }
Real time
For time-dependent problems, this is the time t at the beginning of the current timestep.
Definition: system.h:1595
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
void remove_vector(std::string_view vec_name)
Removes the additional vector vec_name from this system.
Definition: system.C:846
virtual void retrieve_timestep() override
This method retrieves all the stored solutions at the current system.time.
virtual void adjoint_qoi_parameter_sensitivity(const QoISet &qoi_indices, const ParameterVector &parameters, SensitivityData &sensitivities) override
Solves for the derivative of each of the system&#39;s quantities of interest q in qoi[qoi_indices] with r...
Real deltat
For time-dependent problems, this is the amount delta t to advance the solution in time...
Definition: diff_system.h:278
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
IntRange< T > make_range(T beg, T end)
The 2-parameter make_range() helper function returns an IntRange<T> when both input parameters are of...
Definition: int_range.h:134

◆ integrate_qoi_timestep()

void libMesh::UnsteadySolver::integrate_qoi_timestep ( )
overridevirtualinherited

A method to integrate the system::QoI functionals.

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::EulerSolver, libMesh::FirstOrderUnsteadySolver, libMesh::AdaptiveTimeSolver, libMesh::TwostepTimeSolver, and libMesh::Euler2Solver.

Definition at line 325 of file unsteady_solver.C.

326 {
327  libmesh_not_implemented();
328 }

◆ 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 277 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

278  { return _is_adjoint; }
bool _is_adjoint
This boolean tells the TimeSolver whether we are solving a primal or adjoint problem.
Definition: time_solver.h:340

◆ is_steady()

virtual bool libMesh::UnsteadySolver::is_steady ( ) const
inlineoverridevirtualinherited

This is not a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 194 of file unsteady_solver.h.

194 { return false; }

◆ last_completed_timestep_size()

Real libMesh::TimeSolver::last_completed_timestep_size ( )
virtualinherited

Returns system.deltat if fixed timestep solver is used, the complete timestep size (sum of all substeps) if the adaptive time solver is used.

Returns the change in system.time, deltat, for the last timestep which was successfully completed. This only returns the outermost step size in the case of nested time solvers. If no time step has yet been successfully completed, then returns system.deltat.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 160 of file time_solver.C.

References libMesh::TimeSolver::last_deltat.

161 {
162  return last_deltat;
163 }
Real last_deltat
The deltat for the last completed timestep before the current one.
Definition: time_solver.h:332

◆ 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 225 of file time_solver.h.

References libMesh::TimeSolver::_linear_solver.

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

225 { return _linear_solver; }
std::unique_ptr< LinearSolver< Number > > _linear_solver
An implicit linear solver to use for adjoint problems.
Definition: time_solver.h:307

◆ n_objects()

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

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

Definition at line 85 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

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

86  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects
The number of objects.

◆ nonlocal_residual()

bool libMesh::NewmarkSolver::nonlocal_residual ( bool  request_jacobian,
DiffContext context 
)
overridevirtual

This method uses the DifferentiablePhysics' nonlocal_time_derivative() and nonlocal_constraint() to build a full residual for non-local terms.

What combination it uses will depend on theta.

Implements libMesh::TimeSolver.

Definition at line 189 of file newmark_solver.C.

References _general_residual(), libMesh::DifferentiablePhysics::nonlocal_constraint(), libMesh::DifferentiablePhysics::nonlocal_damping_residual(), libMesh::DifferentiablePhysics::nonlocal_mass_residual(), libMesh::DiffContext::nonlocal_reinit(), and libMesh::DifferentiablePhysics::nonlocal_time_derivative().

191 {
192  return this->_general_residual(request_jacobian,
193  context,
199 }
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Adds any nonlocal time derivative contributions (e.g.
Definition: diff_physics.h:214
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Adds any nonlocal constraint contributions (e.g.
Definition: diff_physics.h:233
virtual bool _general_residual(bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit)
This method is the underlying implementation of the public residual methods.
virtual bool nonlocal_damping_residual(bool request_jacobian, DiffContext &)
Subtracts any nonlocal damping vector contributions (e.g.
Definition: diff_physics.h:394
virtual bool nonlocal_mass_residual(bool request_jacobian, DiffContext &c)
Subtracts any nonlocal mass vector contributions (e.g.
Definition: diff_physics.C:57
virtual void nonlocal_reinit(Real)
Gives derived classes the opportunity to reinitialize data needed for nonlocal calculations at a new ...
Definition: diff_context.h:95

◆ old_nonlinear_solution()

Number libMesh::UnsteadySolver::old_nonlinear_solution ( const dof_id_type  global_dof_number) const
inherited
Returns
The old nonlinear solution for the specified global DOF.

Definition at line 337 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::n_dofs(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), _general_residual(), and libMesh::FEMPhysics::eulerian_residual().

339 {
340  libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
341  libmesh_assert_less (global_dof_number, old_local_nonlinear_solution->size());
342 
343  return (*old_local_nonlinear_solution)(global_dof_number);
344 }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
dof_id_type n_dofs() const
Definition: dof_map.h:659
std::shared_ptr< NumericVector< Number > > old_local_nonlinear_solution
Serial vector of _system.get_vector("_old_nonlinear_solution") This is a shared_ptr so that it can be...
const DofMap & get_dof_map() const
Definition: system.h:2293

◆ old_solution_accel()

Number libMesh::SecondOrderUnsteadySolver::old_solution_accel ( const dof_id_type  global_dof_number) const
inherited
Returns
The solution acceleration at the previous time step, \(\ddot{u}_n\), for the specified global DOF.

Definition at line 116 of file second_order_unsteady_solver.C.

References libMesh::SecondOrderUnsteadySolver::_old_local_solution_accel, libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), and libMesh::DofMap::n_dofs().

Referenced by _general_residual(), advance_timestep(), project_initial_accel(), and libMesh::SecondOrderUnsteadySolver::reinit().

118 {
119  libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
120  libmesh_assert_less (global_dof_number, _old_local_solution_accel->size());
121 
122  return (*_old_local_solution_accel)(global_dof_number);
123 }
std::unique_ptr< NumericVector< Number > > _old_local_solution_accel
Serial vector of previous time step acceleration .
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
dof_id_type n_dofs() const
Definition: dof_map.h:659
const DofMap & get_dof_map() const
Definition: system.h:2293

◆ old_solution_rate()

Number libMesh::SecondOrderUnsteadySolver::old_solution_rate ( const dof_id_type  global_dof_number) const
inherited
Returns
The solution rate at the previous time step, \(\dot{u}_n\), for the specified global DOF.

Definition at line 107 of file second_order_unsteady_solver.C.

References libMesh::SecondOrderUnsteadySolver::_old_local_solution_rate, libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), and libMesh::DofMap::n_dofs().

Referenced by _general_residual(), advance_timestep(), libMesh::SecondOrderUnsteadySolver::project_initial_rate(), and libMesh::SecondOrderUnsteadySolver::reinit().

109 {
110  libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
111  libmesh_assert_less (global_dof_number, _old_local_solution_rate->size());
112 
113  return (*_old_local_solution_rate)(global_dof_number);
114 }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
dof_id_type n_dofs() const
Definition: dof_map.h:659
std::unique_ptr< NumericVector< Number > > _old_local_solution_rate
Serial vector of previous time step velocity .
const DofMap & get_dof_map() const
Definition: system.h:2293

◆ print_info()

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

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

Definition at line 81 of file reference_counter.C.

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

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

82 {
84  out_stream << ReferenceCounter::get_info();
85 }
static std::string get_info()
Gets a string containing the reference information.
static bool _enable_print_counter
Flag to control whether reference count information is printed when print_info is called...

◆ project_initial_accel()

void libMesh::NewmarkSolver::project_initial_accel ( FunctionBase< Number > *  f,
FunctionBase< Gradient > *  g = nullptr 
)

Specify non-zero initial acceleration.

Should be called before solve(). This is an alternative to compute_initial_acceleration() if the initial acceleration is actually known. The function value f and its gradient g are user-provided cloneable functors. A gradient g is only required/used for projecting onto finite element spaces with continuous derivatives.

Definition at line 133 of file newmark_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::get_vector(), libMesh::SecondOrderUnsteadySolver::old_solution_accel(), libMesh::System::project_vector(), and set_initial_accel_avail().

135 {
137  _system.get_vector("_old_solution_accel");
138 
140 
141  this->set_initial_accel_avail(true);
142 }
Number old_solution_accel(const dof_id_type global_dof_number) const
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
void set_initial_accel_avail(bool initial_accel_set)
Allow the user to (re)set whether the initial acceleration is available.
void project_vector(NumericVector< Number > &new_vector, FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr, int is_adjoint=-1) const
Projects arbitrary functions onto a vector of degree of freedom values for the current system...
template class LIBMESH_EXPORT NumericVector< Number >
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918

◆ project_initial_rate()

void libMesh::SecondOrderUnsteadySolver::project_initial_rate ( FunctionBase< Number > *  f,
FunctionBase< Gradient > *  g = nullptr 
)
inherited

Specify non-zero initial velocity.

Should be called before solve(). The function value f and its gradient g are user-provided cloneable functors. A gradient g is only required/used for projecting onto finite element spaces with continuous derivatives.

Definition at line 98 of file second_order_unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::get_vector(), libMesh::SecondOrderUnsteadySolver::old_solution_rate(), and libMesh::System::project_vector().

100 {
102  _system.get_vector("_old_solution_rate");
103 
104  _system.project_vector( old_solution_rate, f, g );
105 }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
Number old_solution_rate(const dof_id_type global_dof_number) const
void project_vector(NumericVector< Number > &new_vector, FunctionBase< Number > *f, FunctionBase< Gradient > *g=nullptr, int is_adjoint=-1) const
Projects arbitrary functions onto a vector of degree of freedom values for the current system...
template class LIBMESH_EXPORT NumericVector< Number >
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918

◆ reinit()

void libMesh::SecondOrderUnsteadySolver::reinit ( )
overridevirtualinherited

The reinitialization function.

This method is used to resize internal data vectors after a mesh change.

Reimplemented from libMesh::UnsteadySolver.

Definition at line 60 of file second_order_unsteady_solver.C.

References libMesh::SecondOrderUnsteadySolver::_old_local_solution_accel, libMesh::SecondOrderUnsteadySolver::_old_local_solution_rate, libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::GHOSTED, libMesh::System::n_dofs(), libMesh::System::n_local_dofs(), libMesh::SecondOrderUnsteadySolver::old_solution_accel(), libMesh::SecondOrderUnsteadySolver::old_solution_rate(), libMesh::UnsteadySolver::reinit(), and libMesh::SERIAL.

61 {
63 
64 #ifdef LIBMESH_ENABLE_GHOSTED
67  GHOSTED);
68 
71  GHOSTED);
72 #else
75 #endif
76 
77  // localize the old solutions
79  _system.get_vector("_old_solution_rate");
80 
82  _system.get_vector("_old_solution_accel");
83 
84  old_solution_rate.localize
87 
88  old_solution_accel.localize
91 }
virtual void reinit() override
The reinitialization function.
Number old_solution_accel(const dof_id_type global_dof_number) const
std::unique_ptr< NumericVector< Number > > _old_local_solution_accel
Serial vector of previous time step acceleration .
dof_id_type n_local_dofs() const
Definition: system.C:150
dof_id_type n_dofs() const
Definition: system.C:113
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
Number old_solution_rate(const dof_id_type global_dof_number) const
std::unique_ptr< NumericVector< Number > > _old_local_solution_rate
Serial vector of previous time step velocity .
const DofMap & get_dof_map() const
Definition: system.h:2293
template class LIBMESH_EXPORT NumericVector< Number >
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:511
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918

◆ retrieve_timestep()

void libMesh::SecondOrderUnsteadySolver::retrieve_timestep ( )
overridevirtualinherited

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

Reimplemented from libMesh::UnsteadySolver.

Definition at line 93 of file second_order_unsteady_solver.C.

94 {
95  libmesh_not_implemented();
96 }

◆ set_beta()

void libMesh::NewmarkSolver::set_beta ( Real  beta)
inline

Setter for \( \beta \).

Definition at line 149 of file newmark_solver.h.

References _beta.

150  { _beta = beta; }
Real _beta
The value for the to employ.

◆ set_first_adjoint_step()

void libMesh::UnsteadySolver::set_first_adjoint_step ( bool  first_adjoint_step_setting)
inlineinherited

A setter for the first_adjoint_step boolean.

Needed for nested time solvers.

Definition at line 199 of file unsteady_solver.h.

References libMesh::UnsteadySolver::first_adjoint_step.

200  {
201  first_adjoint_step = first_adjoint_step_setting;
202  }
bool first_adjoint_step
A bool that will be true the first time adjoint_advance_timestep() is called, (when the primal soluti...

◆ set_first_solve()

void libMesh::UnsteadySolver::set_first_solve ( bool  first_solve_setting)
inlineinherited

Definition at line 209 of file unsteady_solver.h.

References libMesh::UnsteadySolver::first_solve.

210  {
211  first_solve = first_solve_setting;
212  }
bool first_solve
A bool that will be true the first time solve() is called, and false thereafter.

◆ set_gamma()

void libMesh::NewmarkSolver::set_gamma ( Real  gamma)
inline

Setter for \( \gamma \).

Note
The method is second order only for \( \gamma = 1/2 \).

Definition at line 157 of file newmark_solver.h.

References _gamma.

158  { _gamma = gamma; }
Real _gamma
The value for to employ.

◆ set_initial_accel_avail()

void libMesh::NewmarkSolver::set_initial_accel_avail ( bool  initial_accel_set)

Allow the user to (re)set whether the initial acceleration is available.

This is not needed if either compute_initial_accel() or project_initial_accel() are called. This is useful is the user is restarting a calculation and the acceleration is available from the restart.

Definition at line 144 of file newmark_solver.C.

References _initial_accel_set.

Referenced by compute_initial_accel(), and project_initial_accel().

145 {
146  _initial_accel_set = initial_accel_set;
147 }
bool _initial_accel_set
This method requires an initial acceleration.

◆ 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 284 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

285  { _is_adjoint = _is_adjoint_value; }
bool _is_adjoint
This boolean tells the TimeSolver whether we are solving a primal or adjoint problem.
Definition: time_solver.h:340

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

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

Referenced by libMesh::AdaptiveTimeSolver::init().

120 {
121  solution_history = _solution_history.clone();
122 }
std::unique_ptr< SolutionHistory > solution_history
A std::unique_ptr to a SolutionHistory object.
Definition: time_solver.h:319

◆ side_residual()

bool libMesh::NewmarkSolver::side_residual ( bool  request_jacobian,
DiffContext context 
)
overridevirtual

This method uses the DifferentiablePhysics' side_time_derivative() and side_constraint() to build a full residual on an element's side.

What combination it uses will depend on theta.

Implements libMesh::TimeSolver.

Definition at line 175 of file newmark_solver.C.

References _general_residual(), libMesh::DiffContext::elem_side_reinit(), libMesh::DifferentiablePhysics::side_constraint(), libMesh::DifferentiablePhysics::side_damping_residual(), libMesh::DifferentiablePhysics::side_mass_residual(), and libMesh::DifferentiablePhysics::side_time_derivative().

177 {
178  return this->_general_residual(request_jacobian,
179  context,
185 }
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Adds the constraint contribution on side of elem to elem_residual.
Definition: diff_physics.h:195
virtual bool side_damping_residual(bool request_jacobian, DiffContext &)
Subtracts a damping vector contribution on side of elem from elem_residual.
Definition: diff_physics.h:378
virtual bool _general_residual(bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit)
This method is the underlying implementation of the public residual methods.
virtual void elem_side_reinit(Real)
Gives derived classes the opportunity to reinitialize data needed for a side integration at a new poi...
Definition: diff_context.h:83
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

◆ solve()

void libMesh::NewmarkSolver::solve ( )
overridevirtual

This method solves for the solution at the next timestep.

Usually we will only need to solve one (non)linear system per timestep, but more complex subclasses may override this.

Reimplemented from libMesh::UnsteadySolver.

Definition at line 149 of file newmark_solver.C.

References _initial_accel_set, and libMesh::UnsteadySolver::solve().

150 {
151  // First, check that the initial accel was set one way or another
152  libmesh_error_msg_if(!_initial_accel_set,
153  "ERROR: Must first set initial acceleration using one of:\n"
154  "NewmarkSolver::compute_initial_accel()\n"
155  "NewmarkSolver::project_initial_accel()\n");
156 
157  // That satisfied, now we can solve
159 }
bool _initial_accel_set
This method requires an initial acceleration.
virtual void solve() override
This method solves for the solution at the next timestep.

◆ system() [1/2]

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

Definition at line 210 of file time_solver.h.

References libMesh::TimeSolver::_system.

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

210 { return _system; }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312

◆ system() [2/2]

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

Definition at line 215 of file time_solver.h.

References libMesh::TimeSolver::_system.

215 { return _system; }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312

◆ time_order()

virtual unsigned int libMesh::SecondOrderUnsteadySolver::time_order ( ) const
inlineoverridevirtualinherited
Returns
The maximum order of time derivatives for which the UnsteadySolver subclass is capable of handling.

For example, EulerSolver will have time_order() = 1 and NewmarkSolver will have time_order() = 2.

Implements libMesh::UnsteadySolver.

Definition at line 60 of file second_order_unsteady_solver.h.

61  { return 2; }

◆ update()

void libMesh::UnsteadySolver::update ( )
inherited

Definition at line 361 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::NumericVector< T >::localize(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.

362 {
363  // Dont forget to localize the old_nonlinear_solution !
364  _system.get_vector("_old_nonlinear_solution").localize
367 }
sys_type & _system
A reference to the system we are solving.
Definition: time_solver.h:312
std::shared_ptr< NumericVector< Number > > old_local_nonlinear_solution
Serial vector of _system.get_vector("_old_nonlinear_solution") This is a shared_ptr so that it can be...
const DofMap & get_dof_map() const
Definition: system.h:2293
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:511
const NumericVector< Number > & get_vector(std::string_view vec_name) const
Definition: system.C:918
virtual void localize(std::vector< T > &v_local) const =0
Creates a copy of the global vector in the local vector v_local.

Member Data Documentation

◆ _beta

Real libMesh::NewmarkSolver::_beta
protected

The value for the \(\beta\) to employ.

Method is unconditionally stable for \( \beta \ge \frac{1}{4} \left( \gamma + \frac{1}{2}\right)^2 \)

Definition at line 167 of file newmark_solver.h.

Referenced by _general_residual(), advance_timestep(), and set_beta().

◆ _counts

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

Actually holds the data.

Definition at line 124 of file reference_counter.h.

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

◆ _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 302 of file time_solver.h.

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

◆ _gamma

Real libMesh::NewmarkSolver::_gamma
protected

The value for \(\gamma\) to employ.

Newmark is 2nd order iff \(\gamma=0.5\).

Definition at line 173 of file newmark_solver.h.

Referenced by _general_residual(), advance_timestep(), error_order(), and set_gamma().

◆ _initial_accel_set

bool libMesh::NewmarkSolver::_initial_accel_set
protected

This method requires an initial acceleration.

So, we force the user to call either compute_initial_accel or project_initial_accel to set the initial acceleration.

Definition at line 187 of file newmark_solver.h.

Referenced by set_initial_accel_avail(), and solve().

◆ _is_accel_solve

bool libMesh::NewmarkSolver::_is_accel_solve
protected

Need to be able to indicate to _general_residual if we are doing an acceleration solve or not.

Definition at line 179 of file newmark_solver.h.

Referenced by _general_residual(), and compute_initial_accel().

◆ _linear_solver

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

An implicit linear solver to use for adjoint problems.

Definition at line 307 of file time_solver.h.

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

◆ _mutex

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.

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

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

◆ _old_local_solution_accel

std::unique_ptr<NumericVector<Number> > libMesh::SecondOrderUnsteadySolver::_old_local_solution_accel
protectedinherited

◆ _old_local_solution_rate

std::unique_ptr<NumericVector<Number> > libMesh::SecondOrderUnsteadySolver::_old_local_solution_rate
protectedinherited

◆ _system

sys_type& libMesh::TimeSolver::_system
protectedinherited

A reference to the system we are solving.

Definition at line 312 of file time_solver.h.

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

◆ first_adjoint_step

bool libMesh::UnsteadySolver::first_adjoint_step
protectedinherited

A bool that will be true the first time adjoint_advance_timestep() is called, (when the primal solution is to be used to set adjoint boundary conditions) and false thereafter.

Definition at line 226 of file unsteady_solver.h.

Referenced by libMesh::AdaptiveTimeSolver::adjoint_advance_timestep(), and libMesh::UnsteadySolver::set_first_adjoint_step().

◆ first_solve

bool libMesh::UnsteadySolver::first_solve
protectedinherited

◆ last_deltat

Real libMesh::TimeSolver::last_deltat
protectedinherited

◆ last_step_deltat

Real libMesh::UnsteadySolver::last_step_deltat
protectedinherited

We will need to move the system.time around to ensure that residuals are built with the right deltat and the right time.

Definition at line 237 of file unsteady_solver.h.

Referenced by libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), and libMesh::EulerSolver::integrate_adjoint_refinement_error_estimate().

◆ next_step_deltat

Real libMesh::UnsteadySolver::next_step_deltat
protectedinherited

◆ old_adjoints

std::vector< std::unique_ptr<NumericVector<Number> > > libMesh::UnsteadySolver::old_adjoints
protectedinherited

A vector of pointers to vectors holding the adjoint solution at the last time step.

Definition at line 231 of file unsteady_solver.h.

Referenced by libMesh::Euler2Solver::integrate_adjoint_refinement_error_estimate(), libMesh::EulerSolver::integrate_adjoint_refinement_error_estimate(), and libMesh::UnsteadySolver::UnsteadySolver().

◆ old_local_nonlinear_solution

std::shared_ptr<NumericVector<Number> > libMesh::UnsteadySolver::old_local_nonlinear_solution
inherited

◆ quiet

bool libMesh::TimeSolver::quiet
inherited

Print extra debugging information if quiet == false.

Definition at line 230 of file time_solver.h.

Referenced by libMesh::TwostepTimeSolver::solve(), libMesh::UnsteadySolver::solve(), and libMesh::EigenTimeSolver::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 259 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

The documentation for this class was generated from the following files: