libMesh::AdaptiveTimeSolver Class Referenceabstract

#include <adaptive_time_solver.h>

Inheritance diagram for libMesh::AdaptiveTimeSolver:

Public Types

typedef FirstOrderUnsteadySolver Parent
 
typedef DifferentiableSystem sys_type
 

Public Member Functions

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

Static Public Member Functions

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

Public Attributes

UniquePtr< UnsteadySolvercore_time_solver
 
SystemNorm component_norm
 
std::vector< float > component_scale
 
Real target_tolerance
 
Real upper_tolerance
 
Real max_deltat
 
Real min_deltat
 
Real max_growth
 
bool global_tolerance
 
UniquePtr< 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 Real calculate_norm (System &, NumericVector< Number > &)
 
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

Real last_deltat
 
bool first_solve
 
bool first_adjoint_step
 
UniquePtr< DiffSolver_diff_solver
 
UniquePtr< LinearSolver< Number > > _linear_solver
 
sys_type_system
 
UniquePtr< SolutionHistorysolution_history
 

Static Protected Attributes

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

Detailed Description

This class wraps another UnsteadySolver derived class, and compares the results of timestepping with deltat and timestepping with 2*deltat to adjust future timestep lengths.

Currently this class only works on fully coupled Systems

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
2007

Definition at line 49 of file adaptive_time_solver.h.

Member Typedef Documentation

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

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

Definition at line 110 of file reference_counter.h.

The parent class

Definition at line 55 of file adaptive_time_solver.h.

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

Definition at line 271 of file time_solver.h.

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

Definitions of argument types for use in refactoring subclasses.

Definition at line 269 of file time_solver.h.

The type of system

Definition at line 64 of file time_solver.h.

Constructor & Destructor Documentation

libMesh::AdaptiveTimeSolver::AdaptiveTimeSolver ( sys_type s)
explicit

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

Definition at line 27 of file adaptive_time_solver.C.

References libMesh::UnsteadySolver::old_local_nonlinear_solution.

30  target_tolerance(1.e-3),
31  upper_tolerance(0.0),
32  max_deltat(0.),
33  min_deltat(0.),
34  max_growth(0.),
35  global_tolerance(true)
36 {
37  // the child class must populate core_time_solver
38  // with whatever actual time solver is to be used
39 
40  // As an UnsteadySolver, we have an old_local_nonlinear_solution, but we're
41  // going to drop it and use our core time solver's instead.
43 }
UniquePtr< UnsteadySolver > core_time_solver
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
libMesh::AdaptiveTimeSolver::~AdaptiveTimeSolver ( )
virtual

Destructor.

Definition at line 47 of file adaptive_time_solver.C.

References libMesh::UnsteadySolver::old_local_nonlinear_solution.

48 {
49  // As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
50  // is being managed by our core_time_solver. Make sure we don't delete
51  // it out from under them, in case the user wants to keep using the core
52  // solver after they're done with us.
54 }
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution

Member Function Documentation

void libMesh::UnsteadySolver::adjoint_advance_timestep ( )
virtualinherited

This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will be done before every UnsteadySolver::adjoint_solve().

Reimplemented 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 }
sys_type & _system
Definition: time_solver.h:256
const DofMap & get_dof_map() const
Definition: system.h:2019
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:788
UniquePtr< SolutionHistory > solution_history
Definition: time_solver.h:263
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:392
virtual void localize(std::vector< T > &v_local) const =0
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
void libMesh::AdaptiveTimeSolver::advance_timestep ( )
virtual

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 87 of file adaptive_time_solver.C.

References libMesh::TimeSolver::_system, libMesh::UnsteadySolver::first_solve, libMesh::System::get_vector(), last_deltat, libMesh::System::solution, and libMesh::System::time.

88 {
89  NumericVector<Number> & old_nonlinear_soln =
90  _system.get_vector("_old_nonlinear_solution");
91  NumericVector<Number> & nonlinear_solution =
92  *(_system.solution);
93  // _system.get_vector("_nonlinear_solution");
94 
95  old_nonlinear_soln = nonlinear_solution;
96 
97  if (!first_solve)
99 }
sys_type & _system
Definition: time_solver.h:256
UniquePtr< NumericVector< Number > > solution
Definition: system.h:1512
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:788
virtual void libMesh::TimeSolver::before_timestep ( )
inlinevirtualinherited

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

Definition at line 166 of file time_solver.h.

166 {}
Real libMesh::AdaptiveTimeSolver::calculate_norm ( System s,
NumericVector< Number > &  v 
)
protectedvirtual

A helper function to calculate error norms

Definition at line 156 of file adaptive_time_solver.C.

References libMesh::System::calculate_norm(), and component_norm.

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

158 {
159  return s.calculate_norm(v, component_norm);
160 }
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::libmesh_assert(), libmesh_nullptr, libMesh::QBase::n_points(), libMesh::DiffContext::n_vars(), libMesh::Variable::type(), and libMesh::System::variable().

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

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 = libmesh_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 }
const class libmesh_nullptr_t libmesh_nullptr
libmesh_assert(j)
unsigned int get_second_order_dot_var(unsigned int var) const
Definition: diff_system.C:311
sys_type & _system
Definition: time_solver.h:256
FEGenericBase< Real > FEBase
bool is_second_order_var(unsigned int var) const
Definition: diff_physics.h:531
UniquePtr< DiffSolver > & libMesh::AdaptiveTimeSolver::diff_solver ( )
virtual

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

Reimplemented from libMesh::TimeSolver.

Definition at line 142 of file adaptive_time_solver.C.

References core_time_solver.

143 {
144  return core_time_solver->diff_solver();
145 }
UniquePtr< UnsteadySolver > core_time_solver
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited
Real libMesh::UnsteadySolver::du ( const SystemNorm norm) const
virtualinherited

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

Note that, 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  UniquePtr<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:256
UniquePtr< NumericVector< Number > > solution
Definition: system.h:1512
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:788
Real calculate_norm(const NumericVector< Number > &v, unsigned int var, FEMNormType norm_type, std::set< unsigned int > *skip_dimensions=libmesh_nullptr) const
Definition: system.C:1393
bool libMesh::AdaptiveTimeSolver::element_residual ( bool  get_jacobian,
DiffContext context 
)
virtual

This method is passed on to the core_time_solver

Implements libMesh::TimeSolver.

Definition at line 112 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::libmesh_assert().

114 {
116 
117  return core_time_solver->element_residual(request_jacobian, context);
118 }
libmesh_assert(j)
UniquePtr< UnsteadySolver > core_time_solver
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.

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

101 {
102  _enable_print_counter = true;
103  return;
104 }
Real libMesh::AdaptiveTimeSolver::error_order ( ) const
virtual

This method is passed on to the core_time_solver

Implements libMesh::UnsteadySolver.

Definition at line 103 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::libmesh_assert().

104 {
106 
107  return core_time_solver->error_order();
108 }
libmesh_assert(j)
UniquePtr< UnsteadySolver > core_time_solver
std::string libMesh::ReferenceCounter::get_info ( )
staticinherited

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

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

48 {
49 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
50 
51  std::ostringstream oss;
52 
53  oss << '\n'
54  << " ---------------------------------------------------------------------------- \n"
55  << "| Reference count information |\n"
56  << " ---------------------------------------------------------------------------- \n";
57 
58  for (Counts::iterator it = _counts.begin();
59  it != _counts.end(); ++it)
60  {
61  const std::string name(it->first);
62  const unsigned int creations = it->second.first;
63  const unsigned int destructions = it->second.second;
64 
65  oss << "| " << name << " reference count information:\n"
66  << "| Creations: " << creations << '\n'
67  << "| Destructions: " << destructions << '\n';
68  }
69 
70  oss << " ---------------------------------------------------------------------------- \n";
71 
72  return oss.str();
73 
74 #else
75 
76  return "";
77 
78 #endif
79 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name)
inlineprotectedinherited

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

Definition at line 160 of file reference_counter.h.

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

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

161 {
162  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
163  std::pair<unsigned int, unsigned int> & p = _counts[name];
164 
165  p.first++;
166 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

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

Definition at line 173 of file reference_counter.h.

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

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

174 {
175  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
176  std::pair<unsigned int, unsigned int> & p = _counts[name];
177 
178  p.second++;
179 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::AdaptiveTimeSolver::init ( )
virtual

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

Reimplemented from libMesh::UnsteadySolver.

Definition at line 58 of file adaptive_time_solver.C.

References core_time_solver, libMesh::libmesh_assert(), and libMesh::UnsteadySolver::old_local_nonlinear_solution.

59 {
61 
62  // We override this because our core_time_solver is the one that
63  // needs to handle new vectors, diff_solver->init(), etc
64  core_time_solver->init();
65 
66  // As an UnsteadySolver, we have an old_local_nonlinear_solution, but it
67  // isn't pointing to the right place - fix it
68  //
69  // This leaves us with two UniquePtrs holding the same pointer - dangerous
70  // for future use. Replace with shared_ptr?
72  UniquePtr<NumericVector<Number> >(core_time_solver->old_local_nonlinear_solution.get());
73 }
libmesh_assert(j)
UniquePtr< UnsteadySolver > core_time_solver
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
void libMesh::UnsteadySolver::init_data ( )
virtualinherited

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 }
virtual void init_data()
Definition: time_solver.C:77
sys_type & _system
Definition: time_solver.h:256
const DofMap & get_dof_map() const
Definition: system.h:2019
dof_id_type n_local_dofs() const
Definition: system.C:180
dof_id_type n_dofs() const
Definition: system.C:143
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:392
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
bool libMesh::TimeSolver::is_adjoint ( ) const
inlineinherited

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

Definition at line 231 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

232  { return _is_adjoint; }
virtual bool libMesh::UnsteadySolver::is_steady ( ) const
inlinevirtualinherited

This is not a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 151 of file unsteady_solver.h.

151 { return false; }
UniquePtr< LinearSolver< Number > > & libMesh::AdaptiveTimeSolver::linear_solver ( )
virtual

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

Reimplemented from libMesh::TimeSolver.

Definition at line 149 of file adaptive_time_solver.C.

References core_time_solver.

150 {
151  return core_time_solver->linear_solver();
152 }
UniquePtr< UnsteadySolver > core_time_solver
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited
bool libMesh::AdaptiveTimeSolver::nonlocal_residual ( bool  get_jacobian,
DiffContext context 
)
virtual

This method is passed on to the core_time_solver

Implements libMesh::TimeSolver.

Definition at line 132 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::libmesh_assert().

134 {
136 
137  return core_time_solver->nonlocal_residual(request_jacobian, context);
138 }
libmesh_assert(j)
UniquePtr< UnsteadySolver > core_time_solver
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 libMesh::EulerSolver::_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 }
dof_id_type n_dofs() const
Definition: dof_map.h:506
sys_type & _system
Definition: time_solver.h:256
const DofMap & get_dof_map() const
Definition: system.h:2019
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
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 libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), and libMesh::FirstOrderUnsteadySolver::time_order().

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 }
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

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

Definition at line 88 of file reference_counter.C.

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

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

89 {
91 }
static std::string get_info()
void libMesh::AdaptiveTimeSolver::reinit ( )
virtual

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

Reimplemented from libMesh::UnsteadySolver.

Definition at line 77 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::libmesh_assert().

78 {
80 
81  // We override this because our core_time_solver is the one that
82  // needs to handle new vectors, diff_solver->reinit(), etc
83  core_time_solver->reinit();
84 }
libmesh_assert(j)
UniquePtr< UnsteadySolver > core_time_solver
void libMesh::UnsteadySolver::retrieve_timestep ( )
virtualinherited

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 }
sys_type & _system
Definition: time_solver.h:256
const DofMap & get_dof_map() const
Definition: system.h:2019
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:788
UniquePtr< SolutionHistory > solution_history
Definition: time_solver.h:263
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:392
virtual void localize(std::vector< T > &v_local) const =0
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
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 238 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

239  { _is_adjoint = _is_adjoint_value; }
void libMesh::TimeSolver::set_solution_history ( const SolutionHistory _solution_history)
inherited

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

Definition at line 97 of file time_solver.C.

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

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

This method is passed on to the core_time_solver

Implements libMesh::TimeSolver.

Definition at line 122 of file adaptive_time_solver.C.

References core_time_solver, and libMesh::libmesh_assert().

124 {
126 
127  return core_time_solver->side_residual(request_jacobian, context);
128 }
libmesh_assert(j)
UniquePtr< UnsteadySolver > core_time_solver
virtual void libMesh::AdaptiveTimeSolver::solve ( )
pure virtual

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.

Implemented in libMesh::TwostepTimeSolver.

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

Definition at line 171 of file time_solver.h.

References libMesh::TimeSolver::_system.

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

171 { return _system; }
sys_type & _system
Definition: time_solver.h:256
sys_type& libMesh::TimeSolver::system ( )
inlineinherited
Returns
a writeable reference to the system we are solving.

Definition at line 176 of file time_solver.h.

References libMesh::TimeSolver::_system.

176 { return _system; }
sys_type & _system
Definition: time_solver.h:256
virtual unsigned int libMesh::FirstOrderUnsteadySolver::time_order ( ) const
inlinevirtualinherited

Returns the maximum order of time derivatives for which the UnsteadySolver subclass is capable of handling. E.g. EulerSolver will have time_order = 1 and NewmarkSolver will have time_order = 2

Implements libMesh::UnsteadySolver.

Definition at line 88 of file first_order_unsteady_solver.h.

References libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), and libMesh::FirstOrderUnsteadySolver::prepare_accel().

89  { return 1; }

Member Data Documentation

ReferenceCounter::Counts libMesh::ReferenceCounter::_counts
staticprotectedinherited
UniquePtr<DiffSolver> libMesh::TimeSolver::_diff_solver
protectedinherited

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

Definition at line 246 of file time_solver.h.

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

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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

Definition at line 134 of file reference_counter.h.

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

UniquePtr<LinearSolver<Number> > libMesh::TimeSolver::_linear_solver
protectedinherited

An implicit linear solver to use for adjoint problems.

Definition at line 251 of file time_solver.h.

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

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 128 of file reference_counter.h.

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

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

Definition at line 123 of file reference_counter.h.

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

sys_type& libMesh::TimeSolver::_system
protectedinherited

A reference to the system we are solving.

Definition at line 256 of file time_solver.h.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::NewmarkSolver::advance_timestep(), advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::NewmarkSolver::compute_initial_accel(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::UnsteadySolver::du(), libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), 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(), 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(), 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().

SystemNorm libMesh::AdaptiveTimeSolver::component_norm

Error calculations are done in this norm, DISCRETE_L2 by default.

Definition at line 119 of file adaptive_time_solver.h.

Referenced by calculate_norm().

std::vector<float> libMesh::AdaptiveTimeSolver::component_scale

If component_norms is non-empty, each variable's contribution to the error of a system will also be scaled by component_scale[var], unless component_scale is empty in which case all variables will be weighted equally.

Definition at line 127 of file adaptive_time_solver.h.

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

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

bool libMesh::UnsteadySolver::first_solve
protectedinherited

A bool that will be true the first time solve() is called, and false thereafter

Definition at line 159 of file unsteady_solver.h.

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

bool libMesh::AdaptiveTimeSolver::global_tolerance

This flag, which is true by default, grows (shrinks) the timestep based on the expected global accuracy of the timestepping scheme. Global in this sense means the cumulative final-time accuracy of the scheme. For example, the backward Euler scheme's truncation error is locally of order 2, so that after N timesteps of size deltat, the result is first-order accurate. If you set this to false, you can grow (shrink) your timestep based on the local accuracy rather than the global accuracy of the core TimeSolver. Note that by setting this value to false you may fail to achieve the predicted convergence in time of the underlying method, however it may be possible to get more fine-grained control over step sizes as well.

Definition at line 197 of file adaptive_time_solver.h.

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

Real libMesh::AdaptiveTimeSolver::last_deltat
protected

We need to store the value of the last deltat used, so that advance_timestep() will increment the system time correctly.

Definition at line 206 of file adaptive_time_solver.h.

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

Real libMesh::AdaptiveTimeSolver::max_deltat

Do not allow the adaptive time solver to select deltat > max_deltat. If you use the default max_deltat=0.0, then deltat is unlimited.

Definition at line 167 of file adaptive_time_solver.h.

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

Real libMesh::AdaptiveTimeSolver::max_growth

Do not allow the adaptive time solver to select a new deltat greater than max_growth times the old deltat. If you use the default max_growth=0.0, then the deltat growth is unlimited.

Definition at line 181 of file adaptive_time_solver.h.

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

Real libMesh::AdaptiveTimeSolver::min_deltat

Do not allow the adaptive time solver to select deltat < min_deltat. The default value is 0.0.

Definition at line 173 of file adaptive_time_solver.h.

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

bool libMesh::TimeSolver::quiet
inherited

Print extra debugging information if quiet == false.

Definition at line 191 of file time_solver.h.

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

unsigned int libMesh::TimeSolver::reduce_deltat_on_diffsolver_failure
inherited

This value (which defaults to zero) is the number of times the TimeSolver is allowed to halve deltat and let the DiffSolver repeat the latest failed solve with a reduced timestep. Note that this has no effect for SteadySolvers. Note that 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 219 of file time_solver.h.

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

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

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

Definition at line 263 of file time_solver.h.

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

Real libMesh::AdaptiveTimeSolver::target_tolerance

This tolerance is the target relative error between an exact time integration and a single time step output, scaled by deltat. integrator, scaled by deltat. If the estimated error exceeds or undershoots the target error tolerance, future timesteps will be run with deltat shrunk or grown to compensate.

The default value is 1.0e-2; obviously users should select their own tolerance.

If a negative target_tolerance is specified, then its absolute value is used to scale the estimated error from the first simulation time step, and this becomes the target tolerance of all future time steps.

Definition at line 144 of file adaptive_time_solver.h.

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

Real libMesh::AdaptiveTimeSolver::upper_tolerance

This tolerance is the maximum relative error between an exact time integration and a single time step output, scaled by deltat. If this error tolerance is exceeded by the estimated error of the current time step, that time step will be repeated with a smaller deltat.

If you use the default upper_tolerance=0.0, then the current time step will not be repeated regardless of estimated error.

If a negative upper_tolerance is specified, then its absolute value is used to scale the estimated error from the first simulation time step, and this becomes the upper tolerance of all future time steps.

Definition at line 161 of file adaptive_time_solver.h.

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


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