libMesh::EulerSolver Class Reference

#include <euler_solver.h>

Inheritance diagram for libMesh::EulerSolver:

Public Types

typedef FirstOrderUnsteadySolver Parent
 
typedef DifferentiableSystem sys_type
 

Public Member Functions

 EulerSolver (sys_type &s)
 
virtual ~EulerSolver ()
 
virtual Real error_order () const override
 
virtual bool element_residual (bool request_jacobian, DiffContext &) override
 
virtual bool side_residual (bool request_jacobian, DiffContext &) override
 
virtual bool nonlocal_residual (bool request_jacobian, DiffContext &) override
 
virtual unsigned int time_order () const override
 
virtual void init () override
 
virtual void init_data () override
 
virtual void reinit () override
 
virtual void solve () override
 
virtual void advance_timestep () override
 
virtual void adjoint_advance_timestep () override
 
virtual void retrieve_timestep () override
 
Number old_nonlinear_solution (const dof_id_type global_dof_number) const
 
virtual Real du (const SystemNorm &norm) const override
 
virtual bool is_steady () const override
 
virtual void before_timestep ()
 
const sys_typesystem () const
 
sys_typesystem ()
 
virtual std::unique_ptr< DiffSolver > & diff_solver ()
 
virtual std::unique_ptr< 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

Real theta
 
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
 
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

virtual bool _general_residual (bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit, bool compute_second_order_eqns)
 
void prepare_accel (DiffContext &context)
 
bool compute_second_order_eqns (bool compute_jacobian, DiffContext &c)
 
void increment_constructor_count (const std::string &name)
 
void increment_destructor_count (const std::string &name)
 

Protected Attributes

bool first_solve
 
bool first_adjoint_step
 
std::unique_ptr< DiffSolver_diff_solver
 
std::unique_ptr< LinearSolver< Number > > _linear_solver
 
sys_type_system
 
std::unique_ptr< 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
 

Detailed Description

This class defines a theta-method Euler (defaulting to Backward Euler with theta = 1.0) solver to handle time integration of DifferentiableSystems.

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
Roy H. Stogner
Date
2006

Definition at line 43 of file euler_solver.h.

Member Typedef Documentation

◆ Counts

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

Data structure to log the information. The log is identified by the class name.

Definition at line 117 of file reference_counter.h.

◆ Parent

The parent class

Definition at line 49 of file euler_solver.h.

◆ ReinitFuncType

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

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

◆ sys_type

The type of system

Definition at line 65 of file time_solver.h.

Constructor & Destructor Documentation

◆ EulerSolver()

libMesh::EulerSolver::EulerSolver ( sys_type s)
explicit

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

Definition at line 27 of file euler_solver.C.

◆ ~EulerSolver()

libMesh::EulerSolver::~EulerSolver ( )
virtual

Destructor.

Definition at line 34 of file euler_solver.C.

35 {
36 }

Member Function Documentation

◆ _general_residual()

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

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

Definition at line 99 of file euler_solver.C.

References libMesh::TimeSolver::_system, libMesh::DenseVector< T >::add(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::DifferentiableSystem::deltat, libMesh::DiffContext::elem_solution_derivative, libMesh::DiffContext::elem_solution_rate_derivative, libMesh::DiffContext::fixed_solution_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_rate(), libMesh::DifferentiableSystem::get_physics(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::FirstOrderUnsteadySolver::prepare_accel(), libMesh::DenseVector< T >::size(), libMesh::DenseVector< T >::swap(), libMesh::DenseMatrix< T >::swap(), theta, and libMesh::System::use_fixed_solution.

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

107 {
108  unsigned int n_dofs = context.get_elem_solution().size();
109 
110  // We might need to save the old jacobian in case one of our physics
111  // terms later is unable to update it analytically.
112  DenseMatrix<Number> old_elem_jacobian(n_dofs, n_dofs);
113  if (request_jacobian)
114  old_elem_jacobian.swap(context.get_elem_jacobian());
115 
116  // Local nonlinear solution at old timestep
117  DenseVector<Number> old_elem_solution(n_dofs);
118  for (unsigned int i=0; i != n_dofs; ++i)
119  old_elem_solution(i) =
120  old_nonlinear_solution(context.get_dof_indices()[i]);
121 
122  // Local time derivative of solution
123  context.get_elem_solution_rate() = context.get_elem_solution();
124  context.get_elem_solution_rate() -= old_elem_solution;
125  context.elem_solution_rate_derivative = 1 / _system.deltat;
126  context.get_elem_solution_rate() *=
127  context.elem_solution_rate_derivative;
128 
129  // If we are asked to compute residuals for second order variables,
130  // we also populate the acceleration part so the user can use that.
132  this->prepare_accel(context);
133 
134  // Local nonlinear solution at time t_theta
135  DenseVector<Number> theta_solution(context.get_elem_solution());
136  theta_solution *= theta;
137  theta_solution.add(1. - theta, old_elem_solution);
138 
139  context.elem_solution_derivative = theta;
140  context.fixed_solution_derivative = theta;
141 
142  // If a fixed solution is requested, we'll use theta_solution
144  context.get_elem_fixed_solution() = theta_solution;
145 
146  // Move theta_->elem_, elem_->theta_
147  context.get_elem_solution().swap(theta_solution);
148 
149  // Move the mesh into place first if necessary, set t = t_{\theta}
150  (context.*reinit_func)(theta);
151 
152  // Get the time derivative at t_theta
153  bool jacobian_computed =
154  (_system.get_physics()->*time_deriv)(request_jacobian, context);
155 
156  jacobian_computed = (_system.get_physics()->*mass)(jacobian_computed, context) &&
157  jacobian_computed;
158 
159  // If we have second-order variables, we need to get damping terms
160  // and the velocity equations
162  {
163  jacobian_computed = (_system.get_physics()->*damping)(jacobian_computed, context) &&
164  jacobian_computed;
165 
166  jacobian_computed = this->compute_second_order_eqns(jacobian_computed, context) &&
167  jacobian_computed;
168  }
169 
170  // Restore the elem position if necessary, set t = t_{n+1}
171  (context.*reinit_func)(1);
172 
173  // Move elem_->elem_, theta_->theta_
174  context.get_elem_solution().swap(theta_solution);
175  context.elem_solution_derivative = 1;
176 
177  // Add the constraint term
178  jacobian_computed = (_system.get_physics()->*constraint)(jacobian_computed, context) &&
179  jacobian_computed;
180 
181  // Add back (or restore) the old jacobian
182  if (request_jacobian)
183  {
184  if (jacobian_computed)
185  context.get_elem_jacobian() += old_elem_jacobian;
186  else
187  context.get_elem_jacobian().swap(old_elem_jacobian);
188  }
189 
190  return jacobian_computed;
191 }
bool compute_second_order_eqns(bool compute_jacobian, DiffContext &c)
Number old_nonlinear_solution(const dof_id_type global_dof_number) const
sys_type & _system
Definition: time_solver.h:258
bool use_fixed_solution
Definition: system.h:1493
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:182

◆ adjoint_advance_timestep()

void libMesh::UnsteadySolver::adjoint_advance_timestep ( )
overridevirtualinherited

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::TimeSolver.

Reimplemented in libMesh::NewmarkSolver.

Definition at line 178 of file unsteady_solver.C.

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

179 {
180  // On the first call of this function, we dont save the adjoint solution or
181  // decrement the time, we just call the retrieve function below
182  if (!first_adjoint_step)
183  {
184  // Call the store function to store the last adjoint before decrementing the time
185  solution_history->store();
186  // Decrement the system time
188  }
189  else
190  {
191  first_adjoint_step = false;
192  }
193 
194  // Retrieve the primal solution vectors at this time using the
195  // solution_history object
196  solution_history->retrieve();
197 
198  // Dont forget to localize the old_nonlinear_solution !
199  _system.get_vector("_old_nonlinear_solution").localize
202 }
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
std::unique_ptr< SolutionHistory > solution_history
Definition: time_solver.h:265
sys_type & _system
Definition: time_solver.h:258
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
const DofMap & get_dof_map() const
Definition: system.h:2049
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:450
virtual void localize(std::vector< T > &v_local) const =0

◆ advance_timestep()

void libMesh::UnsteadySolver::advance_timestep ( )
overridevirtualinherited

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::TimeSolver.

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

Definition at line 152 of file unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, libMesh::UnsteadySolver::first_solve, libMesh::System::get_dof_map(), libMesh::DofMap::get_send_list(), libMesh::System::get_vector(), libMesh::NumericVector< T >::localize(), libMesh::UnsteadySolver::old_local_nonlinear_solution, libMesh::System::solution, libMesh::TimeSolver::solution_history, and libMesh::System::time.

Referenced by libMesh::NewmarkSolver::advance_timestep(), and libMesh::UnsteadySolver::solve().

153 {
154  if (!first_solve)
155  {
156  // Store the solution, does nothing by default
157  // User has to attach appropriate solution_history object for this to
158  // actually store anything anywhere
159  solution_history->store();
160 
162  }
163 
164  NumericVector<Number> & old_nonlinear_soln =
165  _system.get_vector("_old_nonlinear_solution");
166  NumericVector<Number> & nonlinear_solution =
167  *(_system.solution);
168 
169  old_nonlinear_soln = nonlinear_solution;
170 
171  old_nonlinear_soln.localize
174 }
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
std::unique_ptr< SolutionHistory > solution_history
Definition: time_solver.h:265
sys_type & _system
Definition: time_solver.h:258
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
std::unique_ptr< NumericVector< Number > > solution
Definition: system.h:1523
const DofMap & get_dof_map() const
Definition: system.h:2049
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:450
virtual void localize(std::vector< T > &v_local) const =0

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

167 {}

◆ compute_second_order_eqns()

bool libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns ( bool  compute_jacobian,
DiffContext c 
)
protectedinherited

If there are second order variables, then we need to compute their residual equations and corresponding Jacobian. The residual equation will simply be $ \dot{u} - v = 0 $, where $ u $ is the second order variable add by the user and $ v $ is the variable added by the time-solver as the "velocity" variable.

Definition at line 32 of file first_order_unsteady_solver.C.

References libMesh::TimeSolver::_system, libMesh::DiffContext::get_dof_indices(), libMesh::DiffContext::get_elem_jacobian(), libMesh::DiffContext::get_elem_residual(), libMesh::DiffContext::get_elem_solution_derivative(), libMesh::DiffContext::get_elem_solution_rate_derivative(), libMesh::FEMContext::get_element_fe(), libMesh::FEMContext::get_element_qrule(), libMesh::DifferentiableSystem::get_second_order_dot_var(), libMesh::DiffContext::get_system(), libMesh::FEMContext::interior_rate(), libMesh::FEMContext::interior_value(), libMesh::DifferentiablePhysics::is_second_order_var(), libMesh::QBase::n_points(), libMesh::DiffContext::n_vars(), libMesh::Variable::type(), and libMesh::System::variable().

Referenced by _general_residual(), libMesh::Euler2Solver::_general_residual(), element_residual(), libMesh::Euler2Solver::element_residual(), nonlocal_residual(), and libMesh::Euler2Solver::nonlocal_residual().

33 {
34  FEMContext & context = cast_ref<FEMContext &>(c);
35 
36  unsigned int n_qpoints = context.get_element_qrule().n_points();
37 
38  for (unsigned int var = 0; var != context.n_vars(); ++var)
39  {
40  if (!this->_system.is_second_order_var(var))
41  continue;
42 
43  unsigned int dot_var = this->_system.get_second_order_dot_var(var);
44 
45  // We're assuming that the FE space for var and dot_var are the same
46  libmesh_assert( context.get_system().variable(var).type() ==
47  context.get_system().variable(dot_var).type() );
48 
49  FEBase * elem_fe = nullptr;
50  context.get_element_fe( var, elem_fe );
51 
52  const std::vector<Real> & JxW = elem_fe->get_JxW();
53 
54  const std::vector<std::vector<Real>> & phi = elem_fe->get_phi();
55 
56  const unsigned int n_dofs = cast_int<unsigned int>
57  (context.get_dof_indices(dot_var).size());
58 
59  DenseSubVector<Number> & Fu = context.get_elem_residual(var);
60  DenseSubMatrix<Number> & Kuu = context.get_elem_jacobian( var, var );
61  DenseSubMatrix<Number> & Kuv = context.get_elem_jacobian( var, dot_var );
62 
63  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
64  {
65  Number udot, v;
66  context.interior_rate(var, qp, udot);
67  context.interior_value(dot_var, qp, v);
68 
69  for (unsigned int i = 0; i < n_dofs; i++)
70  {
71  Fu(i) += JxW[qp]*(udot-v)*phi[i][qp];
72 
73  if (compute_jacobian)
74  {
75  Number rate_factor = JxW[qp]*context.get_elem_solution_rate_derivative()*phi[i][qp];
76  Number soln_factor = JxW[qp]*context.get_elem_solution_derivative()*phi[i][qp];
77 
78  Kuu(i,i) += rate_factor*phi[i][qp];
79  Kuv(i,i) -= soln_factor*phi[i][qp];
80 
81  for (unsigned int j = i+1; j < n_dofs; j++)
82  {
83  Kuu(i,j) += rate_factor*phi[j][qp];
84  Kuu(j,i) += rate_factor*phi[j][qp];
85 
86  Kuv(i,j) -= soln_factor*phi[j][qp];
87  Kuv(j,i) -= soln_factor*phi[j][qp];
88  }
89  }
90  }
91  }
92  }
93 
94  return compute_jacobian;
95 }
bool is_second_order_var(unsigned int var) const
Definition: diff_physics.h:533
sys_type & _system
Definition: time_solver.h:258
unsigned int get_second_order_dot_var(unsigned int var) const
Definition: diff_system.C:306
FEGenericBase< Real > FEBase

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

182 { return _diff_solver; }
std::unique_ptr< DiffSolver > _diff_solver
Definition: time_solver.h:248

◆ disable_print_counter_info()

void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

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

◆ 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 227 of file unsteady_solver.C.

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

228 {
229 
230  std::unique_ptr<NumericVector<Number>> solution_copy =
231  _system.solution->clone();
232 
233  solution_copy->add(-1., _system.get_vector("_old_nonlinear_solution"));
234 
235  solution_copy->close();
236 
237  return _system.calculate_norm(*solution_copy, norm);
238 }
sys_type & _system
Definition: time_solver.h:258
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
std::unique_ptr< NumericVector< Number > > solution
Definition: system.h:1523
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=nullptr) const
Definition: system.C:1378

◆ element_residual()

bool libMesh::EulerSolver::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 50 of file euler_solver.C.

References libMesh::DifferentiablePhysics::_eulerian_time_deriv(), _general_residual(), libMesh::TimeSolver::_system, libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::DifferentiablePhysics::damping_residual(), libMesh::DiffContext::elem_reinit(), libMesh::DifferentiablePhysics::element_constraint(), libMesh::DifferentiableSystem::get_physics(), libMesh::DifferentiablePhysics::have_second_order_vars(), and libMesh::DifferentiablePhysics::mass_residual().

52 {
54 
55  return this->_general_residual(request_jacobian,
56  context,
63 }
bool compute_second_order_eqns(bool compute_jacobian, DiffContext &c)
virtual bool damping_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:374
virtual bool _general_residual(bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit, bool compute_second_order_eqns)
Definition: euler_solver.C:99
virtual void elem_reinit(Real)
Definition: diff_context.h:76
bool _eulerian_time_deriv(bool request_jacobian, DiffContext &)
virtual bool mass_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:318
sys_type & _system
Definition: time_solver.h:258
const DifferentiablePhysics * get_physics() const
Definition: diff_system.h:182
virtual bool element_constraint(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:142

◆ enable_print_counter_info()

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

◆ error_order()

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

Error convergence order: 2 for Crank-Nicolson, 1 otherwise

Implements libMesh::UnsteadySolver.

Definition at line 40 of file euler_solver.C.

References theta.

41 {
42  if (theta == 0.5)
43  return 2.;
44  return 1.;
45 }

◆ 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 & pr : _counts)
59  {
60  const std::string name(pr.first);
61  const unsigned int creations = pr.second.first;
62  const unsigned int destructions = pr.second.second;
63 
64  oss << "| " << name << " reference count information:\n"
65  << "| Creations: " << creations << '\n'
66  << "| Destructions: " << destructions << '\n';
67  }
68 
69  oss << " ---------------------------------------------------------------------------- \n";
70 
71  return oss.str();
72 
73 #else
74 
75  return "";
76 
77 #endif
78 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42

◆ increment_constructor_count()

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

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

Definition at line 181 of file reference_counter.h.

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

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

182 {
183  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
184  std::pair<unsigned int, unsigned int> & p = _counts[name];
185 
186  p.first++;
187 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ increment_destructor_count()

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

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

Definition at line 194 of file reference_counter.h.

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

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

195 {
196  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
197  std::pair<unsigned int, unsigned int> & p = _counts[name];
198 
199  p.second++;
200 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
spin_mutex spin_mtx
Definition: threads.C:29

◆ init()

void libMesh::UnsteadySolver::init ( )
overridevirtualinherited

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

Reimplemented from libMesh::TimeSolver.

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

Definition at line 46 of file unsteady_solver.C.

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

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

47 {
49 
50  _system.add_vector("_old_nonlinear_solution");
51 }
sys_type & _system
Definition: time_solver.h:258
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Definition: system.C:661
virtual void init()
Definition: time_solver.C:64

◆ init_data()

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

Reimplemented in libMesh::SecondOrderUnsteadySolver.

Definition at line 55 of file unsteady_solver.C.

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

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

56 {
58 
59 #ifdef LIBMESH_ENABLE_GHOSTED
62  GHOSTED);
63 #else
65 #endif
66 }
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
virtual void init_data()
Definition: time_solver.C:77
dof_id_type n_local_dofs() const
Definition: system.C:187
dof_id_type n_dofs() const
Definition: system.C:150
sys_type & _system
Definition: time_solver.h:258
const DofMap & get_dof_map() const
Definition: system.h:2049
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:450

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

References libMesh::TimeSolver::_is_adjoint.

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

234  { return _is_adjoint; }

◆ is_steady()

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

This is not a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 154 of file unsteady_solver.h.

154 { return false; }

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

References libMesh::TimeSolver::_linear_solver.

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

187 { return _linear_solver; }
std::unique_ptr< LinearSolver< Number > > _linear_solver
Definition: time_solver.h:253

◆ n_objects()

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

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

Definition at line 83 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

84  { return _n_objects; }
static Threads::atomic< unsigned int > _n_objects

◆ nonlocal_residual()

bool libMesh::EulerSolver::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 82 of file euler_solver.C.

References _general_residual(), libMesh::TimeSolver::_system, libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::DifferentiableSystem::have_second_order_scalar_vars(), libMesh::DifferentiablePhysics::nonlocal_constraint(), libMesh::DifferentiablePhysics::nonlocal_damping_residual(), libMesh::DifferentiablePhysics::nonlocal_mass_residual(), libMesh::DiffContext::nonlocal_reinit(), and libMesh::DifferentiablePhysics::nonlocal_time_derivative().

84 {
86 
87  return this->_general_residual(request_jacobian,
88  context,
95 }
bool compute_second_order_eqns(bool compute_jacobian, DiffContext &c)
virtual bool nonlocal_time_derivative(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:209
virtual bool _general_residual(bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit, bool compute_second_order_eqns)
Definition: euler_solver.C:99
virtual bool nonlocal_constraint(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:227
sys_type & _system
Definition: time_solver.h:258
bool have_second_order_scalar_vars() const
Definition: diff_system.C:335
virtual bool nonlocal_damping_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:406
virtual bool nonlocal_mass_residual(bool request_jacobian, DiffContext &c)
virtual void nonlocal_reinit(Real)
Definition: diff_context.h:94

◆ 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 216 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 _general_residual(), libMesh::Euler2Solver::_general_residual(), and libMesh::NewmarkSolver::_general_residual().

218 {
219  libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
220  libmesh_assert_less (global_dof_number, old_local_nonlinear_solution->size());
221 
222  return (*old_local_nonlinear_solution)(global_dof_number);
223 }
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
sys_type & _system
Definition: time_solver.h:258
dof_id_type n_dofs() const
Definition: dof_map.h:574
const DofMap & get_dof_map() const
Definition: system.h:2049

◆ prepare_accel()

void libMesh::FirstOrderUnsteadySolver::prepare_accel ( DiffContext context)
protectedinherited

If there are second order variables in the system, then we also prepare the accel for those variables so the user can treat them as such.

Definition at line 25 of file first_order_unsteady_solver.C.

References libMesh::DiffContext::elem_solution_accel_derivative, libMesh::DiffContext::get_elem_solution_accel(), libMesh::DiffContext::get_elem_solution_rate(), and libMesh::DiffContext::get_elem_solution_rate_derivative().

Referenced by _general_residual(), and libMesh::Euler2Solver::_general_residual().

26 {
27  context.get_elem_solution_accel() = context.get_elem_solution_rate();
28 
29  context.elem_solution_accel_derivative = context.get_elem_solution_rate_derivative();
30 }

◆ print_info()

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

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

Definition at line 87 of file reference_counter.C.

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

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

◆ reinit()

void libMesh::UnsteadySolver::reinit ( )
overridevirtualinherited

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

Reimplemented from libMesh::TimeSolver.

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

Definition at line 70 of file unsteady_solver.C.

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

Referenced by libMesh::SecondOrderUnsteadySolver::reinit().

71 {
73 
74 #ifdef LIBMESH_ENABLE_GHOSTED
77  GHOSTED);
78 #else
80 #endif
81 
82  // localize the old solution
83  NumericVector<Number> & old_nonlinear_soln =
84  _system.get_vector("_old_nonlinear_solution");
85 
86  old_nonlinear_soln.localize
89 }
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
virtual void reinit()
Definition: time_solver.C:48
dof_id_type n_local_dofs() const
Definition: system.C:187
dof_id_type n_dofs() const
Definition: system.C:150
sys_type & _system
Definition: time_solver.h:258
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
const DofMap & get_dof_map() const
Definition: system.h:2049
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:450
virtual void localize(std::vector< T > &v_local) const =0

◆ retrieve_timestep()

void libMesh::UnsteadySolver::retrieve_timestep ( )
overridevirtualinherited

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

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::SecondOrderUnsteadySolver.

Definition at line 204 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(), libMesh::UnsteadySolver::old_local_nonlinear_solution, and libMesh::TimeSolver::solution_history.

205 {
206  // Retrieve all the stored vectors at the current time
207  solution_history->retrieve();
208 
209  // Dont forget to localize the old_nonlinear_solution !
210  _system.get_vector("_old_nonlinear_solution").localize
213 }
std::unique_ptr< NumericVector< Number > > old_local_nonlinear_solution
std::unique_ptr< SolutionHistory > solution_history
Definition: time_solver.h:265
sys_type & _system
Definition: time_solver.h:258
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:774
const DofMap & get_dof_map() const
Definition: system.h:2049
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:450
virtual void localize(std::vector< T > &v_local) const =0

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

References libMesh::TimeSolver::_is_adjoint.

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

241  { _is_adjoint = _is_adjoint_value; }

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

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

98 {
99  solution_history = _solution_history.clone();
100 }
std::unique_ptr< SolutionHistory > solution_history
Definition: time_solver.h:265

◆ side_residual()

bool libMesh::EulerSolver::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 67 of file euler_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().

69 {
70  return this->_general_residual(request_jacobian,
71  context,
77  false);
78 }
virtual bool side_constraint(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:191
virtual bool side_damping_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:391
virtual bool _general_residual(bool request_jacobian, DiffContext &, ResFuncType mass, ResFuncType damping, ResFuncType time_deriv, ResFuncType constraint, ReinitFuncType reinit, bool compute_second_order_eqns)
Definition: euler_solver.C:99
virtual void elem_side_reinit(Real)
Definition: diff_context.h:82
virtual bool side_mass_residual(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:335
virtual bool side_time_derivative(bool request_jacobian, DiffContext &)
Definition: diff_physics.h:171

◆ solve()

void libMesh::UnsteadySolver::solve ( )
overridevirtualinherited

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::TimeSolver.

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

Definition at line 93 of file unsteady_solver.C.

References libMesh::TimeSolver::_diff_solver, libMesh::TimeSolver::_system, libMesh::UnsteadySolver::advance_timestep(), libMesh::DifferentiableSystem::deltat, libMesh::DiffSolver::DIVERGED_BACKTRACKING_FAILURE, libMesh::DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS, libMesh::UnsteadySolver::first_solve, libMesh::out, libMesh::TimeSolver::quiet, and libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure.

Referenced by libMesh::NewmarkSolver::solve().

94 {
95  if (first_solve)
96  {
98  first_solve = false;
99  }
100 
101  unsigned int solve_result = _diff_solver->solve();
102 
103  // If we requested the UnsteadySolver to attempt reducing dt after a
104  // failed DiffSolver solve, check the results of the solve now.
106  {
107  bool backtracking_failed =
109 
110  bool max_iterations =
112 
113  if (backtracking_failed || max_iterations)
114  {
115  // Cut timestep in half
116  for (unsigned int nr=0; nr<reduce_deltat_on_diffsolver_failure; ++nr)
117  {
118  _system.deltat *= 0.5;
119  libMesh::out << "Newton backtracking failed. Trying with smaller timestep, dt="
120  << _system.deltat << std::endl;
121 
122  solve_result = _diff_solver->solve();
123 
124  // Check solve results with reduced timestep
125  bool backtracking_still_failed =
127 
128  bool backtracking_max_iterations =
130 
131  if (!backtracking_still_failed && !backtracking_max_iterations)
132  {
133  if (!quiet)
134  libMesh::out << "Reduced dt solve succeeded." << std::endl;
135  return;
136  }
137  }
138 
139  // If we made it here, we still couldn't converge the solve after
140  // reducing deltat
141  libMesh::out << "DiffSolver::solve() did not succeed after "
143  << " attempts." << std::endl;
144  libmesh_convergence_failure();
145 
146  } // end if (backtracking_failed || max_iterations)
147  } // end if (reduce_deltat_on_diffsolver_failure)
148 }
std::unique_ptr< DiffSolver > _diff_solver
Definition: time_solver.h:248
sys_type & _system
Definition: time_solver.h:258
unsigned int reduce_deltat_on_diffsolver_failure
Definition: time_solver.h:221
virtual void advance_timestep() override
OStreamProxy out(std::cout)

◆ system() [1/2]

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

Definition at line 172 of file time_solver.h.

References libMesh::TimeSolver::_system.

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

172 { return _system; }
sys_type & _system
Definition: time_solver.h:258

◆ system() [2/2]

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

Definition at line 177 of file time_solver.h.

References libMesh::TimeSolver::_system.

177 { return _system; }
sys_type & _system
Definition: time_solver.h:258

◆ time_order()

virtual unsigned int libMesh::FirstOrderUnsteadySolver::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 90 of file first_order_unsteady_solver.h.

91  { return 1; }

Member Data Documentation

◆ _counts

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

◆ _diff_solver

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

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

Definition at line 248 of file time_solver.h.

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

◆ _enable_print_counter

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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

Definition at line 141 of file reference_counter.h.

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

◆ _linear_solver

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

An implicit linear solver to use for adjoint problems.

Definition at line 253 of file time_solver.h.

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

◆ _mutex

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 135 of file reference_counter.h.

◆ _n_objects

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

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 130 of file reference_counter.h.

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

◆ _system

sys_type& libMesh::TimeSolver::_system
protectedinherited

A reference to the system we are solving.

Definition at line 258 of file time_solver.h.

Referenced by _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(), element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::EigenTimeSolver::element_residual(), libMesh::SecondOrderUnsteadySolver::init(), libMesh::UnsteadySolver::init(), libMesh::TimeSolver::init(), libMesh::EigenTimeSolver::init(), libMesh::SecondOrderUnsteadySolver::init_data(), libMesh::UnsteadySolver::init_data(), libMesh::TimeSolver::init_data(), 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(), libMesh::NewmarkSolver::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(), and libMesh::TimeSolver::system().

◆ 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 168 of file unsteady_solver.h.

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

◆ first_solve

bool libMesh::UnsteadySolver::first_solve
protectedinherited

◆ old_local_nonlinear_solution

◆ quiet

bool libMesh::TimeSolver::quiet
inherited

Print extra debugging information if quiet == false.

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

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

◆ solution_history

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

A std::unique_ptr to a SolutionHistory object. Default is NoSolutionHistory, which the user can override by declaring a different kind of SolutionHistory in the application

Definition at line 265 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().

◆ theta

Real libMesh::EulerSolver::theta

The value for the theta method to employ: 1.0 corresponds to backwards Euler, 0.0 corresponds to forwards Euler, 0.5 corresponds to Crank-Nicolson.

Definition at line 100 of file euler_solver.h.

Referenced by _general_residual(), and error_order().


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