libMesh::SecondOrderUnsteadySolver Class Referenceabstract

#include <second_order_unsteady_solver.h>

Inheritance diagram for libMesh::SecondOrderUnsteadySolver:

Public Types

typedef DifferentiableSystem sys_type
 

Public Member Functions

 SecondOrderUnsteadySolver (sys_type &s)
 
virtual ~SecondOrderUnsteadySolver ()
 
virtual unsigned int time_order () const libmesh_override
 
virtual void init () libmesh_override
 
virtual void init_data () libmesh_override
 
virtual void reinit () libmesh_override
 
virtual void retrieve_timestep () libmesh_override
 
void project_initial_rate (FunctionBase< Number > *f, FunctionBase< Gradient > *g=libmesh_nullptr)
 
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 solve () libmesh_override
 
virtual void advance_timestep () libmesh_override
 
virtual void adjoint_advance_timestep () libmesh_override
 
virtual Real error_order () const =0
 
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 bool element_residual (bool request_jacobian, DiffContext &)=0
 
virtual bool side_residual (bool request_jacobian, DiffContext &)=0
 
virtual bool nonlocal_residual (bool request_jacobian, DiffContext &)=0
 
virtual void before_timestep ()
 
const sys_typesystem () const
 
sys_typesystem ()
 
virtual UniquePtr< DiffSolver > & diff_solver ()
 
virtual UniquePtr< LinearSolver< Number > > & linear_solver ()
 
void set_solution_history (const SolutionHistory &_solution_history)
 
bool is_adjoint () const
 
void set_is_adjoint (bool _is_adjoint_value)
 

Static Public Member Functions

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

Public Attributes

UniquePtr< 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

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

Protected Attributes

UniquePtr< NumericVector< Number > > _old_local_solution_rate
 
UniquePtr< NumericVector< Number > > _old_local_solution_accel
 
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

Generic class from which second order UnsteadySolvers should subclass.

Subclasses of this class are meant to solve problems of the form

\[ M(u)\ddot{u} + C(u)\dot{u} + F(u) = 0 \]

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 38 of file second_order_unsteady_solver.h.

Member Typedef Documentation

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

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

Definition at line 119 of file reference_counter.h.

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

Definition at line 272 of file time_solver.h.

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

Definitions of argument types for use in refactoring subclasses.

Definition at line 270 of file time_solver.h.

The type of system

Definition at line 64 of file time_solver.h.

Constructor & Destructor Documentation

libMesh::SecondOrderUnsteadySolver::SecondOrderUnsteadySolver ( sys_type s)
explicit

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

Definition at line 25 of file second_order_unsteady_solver.C.

26  : UnsteadySolver(s),
29 {}
static UniquePtr< NumericVector< Number > > build(const Parallel::Communicator &comm, const SolverPackage solver_package=libMesh::default_solver_package())
UniquePtr< NumericVector< Number > > _old_local_solution_rate
UniquePtr< NumericVector< Number > > _old_local_solution_accel
libMesh::SecondOrderUnsteadySolver::~SecondOrderUnsteadySolver ( )
virtual

Destructor.

Definition at line 31 of file second_order_unsteady_solver.C.

32 {}

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:257
const DofMap & get_dof_map() const
Definition: system.h:2028
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:792
UniquePtr< SolutionHistory > solution_history
Definition: time_solver.h:264
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:394
virtual void localize(std::vector< T > &v_local) const =0
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
void libMesh::UnsteadySolver::advance_timestep ( )
virtualinherited

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 }
sys_type & _system
Definition: time_solver.h:257
const DofMap & get_dof_map() const
Definition: system.h:2028
UniquePtr< NumericVector< Number > > solution
Definition: system.h:1521
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:792
UniquePtr< SolutionHistory > solution_history
Definition: time_solver.h:264
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:394
virtual void localize(std::vector< T > &v_local) const =0
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
virtual void libMesh::TimeSolver::before_timestep ( )
inlinevirtualinherited

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

Definition at line 166 of file time_solver.h.

166 {}
virtual UniquePtr<DiffSolver>& libMesh::TimeSolver::diff_solver ( )
inlinevirtualinherited

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 181 of file time_solver.h.

References libMesh::TimeSolver::_diff_solver.

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

181 { return _diff_solver; }
UniquePtr< DiffSolver > _diff_solver
Definition: time_solver.h:247
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited
Real libMesh::UnsteadySolver::du ( const SystemNorm norm) const
virtualinherited

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  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:257
UniquePtr< NumericVector< Number > > solution
Definition: system.h:1521
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:792
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:1397
virtual bool libMesh::TimeSolver::element_residual ( bool  request_jacobian,
DiffContext  
)
pure virtualinherited

This method uses the DifferentiablePhysics element_time_derivative(), element_constraint(), and mass_residual() to build a full residual on an element. What combination

it uses will depend on the type of solver. See the subclasses for more details.

Implemented in libMesh::EigenTimeSolver, libMesh::NewmarkSolver, libMesh::AdaptiveTimeSolver, libMesh::SteadySolver, libMesh::Euler2Solver, and libMesh::EulerSolver.

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

void libMesh::ReferenceCounter::enable_print_counter_info ( )
staticinherited

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

Definition at line 101 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

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

102 {
103  _enable_print_counter = true;
104  return;
105 }
virtual Real libMesh::UnsteadySolver::error_order ( ) const
pure virtualinherited

This method should return the expected convergence order of the (non-local) error of the time discretization scheme - e.g. 2 for the O(deltat^2) Crank-Nicholson, or 1 for the O(deltat) Backward Euler.

Useful for adaptive timestepping schemes.

Implemented in libMesh::NewmarkSolver, libMesh::AdaptiveTimeSolver, libMesh::Euler2Solver, and libMesh::EulerSolver.

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

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

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

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

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

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

Definition at line 185 of file reference_counter.h.

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

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

186 {
187  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
188  std::pair<unsigned int, unsigned int> & p = _counts[name];
189 
190  p.first++;
191 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name)
inlineprotectedinherited

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

Definition at line 198 of file reference_counter.h.

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

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

199 {
200  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
201  std::pair<unsigned int, unsigned int> & p = _counts[name];
202 
203  p.second++;
204 }
std::string name(const ElemQuality q)
Definition: elem_quality.C:39
spin_mutex spin_mtx
Definition: threads.C:29
void libMesh::SecondOrderUnsteadySolver::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 34 of file second_order_unsteady_solver.C.

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

Referenced by time_order().

35 {
37 
38  _system.add_vector("_old_solution_rate");
39  _system.add_vector("_old_solution_accel");
40 }
virtual void init() libmesh_override
sys_type & _system
Definition: time_solver.h:257
NumericVector< Number > & add_vector(const std::string &vec_name, const bool projections=true, const ParallelType type=PARALLEL)
Definition: system.C:679
void libMesh::SecondOrderUnsteadySolver::init_data ( )
virtual

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 _old_local_solution_accel, _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.

Referenced by time_order().

43 {
45 
46 #ifdef LIBMESH_ENABLE_GHOSTED
49  GHOSTED);
50 
53  GHOSTED);
54 #else
57 #endif
58 }
sys_type & _system
Definition: time_solver.h:257
const DofMap & get_dof_map() const
Definition: system.h:2028
UniquePtr< NumericVector< Number > > _old_local_solution_rate
dof_id_type n_local_dofs() const
Definition: system.C:183
UniquePtr< NumericVector< Number > > _old_local_solution_accel
dof_id_type n_dofs() const
Definition: system.C:146
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:394
virtual void init_data() libmesh_override
bool libMesh::TimeSolver::is_adjoint ( ) const
inlineinherited

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

Definition at line 232 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

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

This is not a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 153 of file unsteady_solver.h.

153 { return false; }
virtual UniquePtr<LinearSolver<Number> >& libMesh::TimeSolver::linear_solver ( )
inlinevirtualinherited

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 186 of file time_solver.h.

References libMesh::TimeSolver::_linear_solver.

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

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

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

What combination it uses will depend on the type of solver. See the subclasses for more details.

Implemented in libMesh::NewmarkSolver, libMesh::EigenTimeSolver, libMesh::SteadySolver, libMesh::AdaptiveTimeSolver, libMesh::Euler2Solver, and libMesh::EulerSolver.

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

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:510
sys_type & _system
Definition: time_solver.h:257
const DofMap & get_dof_map() const
Definition: system.h:2028
UniquePtr< NumericVector< Number > > old_local_nonlinear_solution
Number libMesh::SecondOrderUnsteadySolver::old_solution_accel ( const dof_id_type  global_dof_number) const
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 _old_local_solution_accel, libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), and libMesh::DofMap::n_dofs().

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

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 }
dof_id_type n_dofs() const
Definition: dof_map.h:510
sys_type & _system
Definition: time_solver.h:257
const DofMap & get_dof_map() const
Definition: system.h:2028
UniquePtr< NumericVector< Number > > _old_local_solution_accel
Number libMesh::SecondOrderUnsteadySolver::old_solution_rate ( const dof_id_type  global_dof_number) const
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 _old_local_solution_rate, libMesh::TimeSolver::_system, libMesh::System::get_dof_map(), and libMesh::DofMap::n_dofs().

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

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 }
dof_id_type n_dofs() const
Definition: dof_map.h:510
sys_type & _system
Definition: time_solver.h:257
const DofMap & get_dof_map() const
Definition: system.h:2028
UniquePtr< NumericVector< Number > > _old_local_solution_rate
void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out)
staticinherited

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

Definition at line 88 of file reference_counter.C.

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

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

89 {
91  out_stream << ReferenceCounter::get_info();
92 }
static std::string get_info()
void libMesh::SecondOrderUnsteadySolver::project_initial_rate ( FunctionBase< Number > *  f,
FunctionBase< Gradient > *  g = libmesh_nullptr 
)

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(), old_solution_rate(), and libMesh::System::project_vector().

Referenced by time_order().

100 {
101  NumericVector<Number> & old_solution_rate =
102  _system.get_vector("_old_solution_rate");
103 
104  _system.project_vector( old_solution_rate, f, g );
105 }
sys_type & _system
Definition: time_solver.h:257
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:792
void project_vector(NumericVector< Number > &new_vector, FunctionBase< Number > *f, FunctionBase< Gradient > *g=libmesh_nullptr, int is_adjoint=-1) const
Number old_solution_rate(const dof_id_type global_dof_number) const
void libMesh::SecondOrderUnsteadySolver::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 60 of file second_order_unsteady_solver.C.

References _old_local_solution_accel, _old_local_solution_rate, 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(), old_solution_accel(), old_solution_rate(), libMesh::UnsteadySolver::reinit(), and libMesh::SERIAL.

Referenced by time_order().

61 {
63 
64 #ifdef LIBMESH_ENABLE_GHOSTED
67  GHOSTED);
68 
71  GHOSTED);
72 #else
75 #endif
76 
77  // localize the old solutions
78  NumericVector<Number> & old_solution_rate =
79  _system.get_vector("_old_solution_rate");
80 
81  NumericVector<Number> & old_solution_accel =
82  _system.get_vector("_old_solution_accel");
83 
84  old_solution_rate.localize
87 
88  old_solution_accel.localize
91 }
Number old_solution_accel(const dof_id_type global_dof_number) const
sys_type & _system
Definition: time_solver.h:257
const DofMap & get_dof_map() const
Definition: system.h:2028
const NumericVector< Number > & get_vector(const std::string &vec_name) const
Definition: system.C:792
UniquePtr< NumericVector< Number > > _old_local_solution_rate
virtual void reinit() libmesh_override
dof_id_type n_local_dofs() const
Definition: system.C:183
UniquePtr< NumericVector< Number > > _old_local_solution_accel
dof_id_type n_dofs() const
Definition: system.C:146
Number old_solution_rate(const dof_id_type global_dof_number) const
const std::vector< dof_id_type > & get_send_list() const
Definition: dof_map.h:394
virtual void localize(std::vector< T > &v_local) const =0
void libMesh::SecondOrderUnsteadySolver::retrieve_timestep ( )
virtual

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.

Referenced by time_order().

94 {
95  libmesh_not_implemented();
96 }
void libMesh::TimeSolver::set_is_adjoint ( bool  _is_adjoint_value)
inlineinherited

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

Definition at line 239 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

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

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

Definition at line 97 of file time_solver.C.

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

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

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

What combination it uses will depend on the type of solver. See the subclasses for more details.

Implemented in libMesh::NewmarkSolver, libMesh::EigenTimeSolver, libMesh::AdaptiveTimeSolver, libMesh::SteadySolver, libMesh::Euler2Solver, and libMesh::EulerSolver.

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

void libMesh::UnsteadySolver::solve ( )
virtualinherited

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 "
142  << reduce_deltat_on_diffsolver_failure
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 }
virtual void advance_timestep() libmesh_override
sys_type & _system
Definition: time_solver.h:257
unsigned int reduce_deltat_on_diffsolver_failure
Definition: time_solver.h:220
UniquePtr< DiffSolver > _diff_solver
Definition: time_solver.h:247
OStreamProxy out(std::cout)
const sys_type& libMesh::TimeSolver::system ( ) const
inlineinherited
Returns
A constant reference to the system we are solving.

Definition at line 171 of file time_solver.h.

References libMesh::TimeSolver::_system.

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

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

Definition at line 176 of file time_solver.h.

References libMesh::TimeSolver::_system.

176 { return _system; }
sys_type & _system
Definition: time_solver.h:257
virtual unsigned int libMesh::SecondOrderUnsteadySolver::time_order ( ) const
inlinevirtual
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 53 of file second_order_unsteady_solver.h.

References init(), init_data(), libmesh_nullptr, old_solution_accel(), old_solution_rate(), project_initial_rate(), reinit(), and retrieve_timestep().

54  { return 2; }

Member Data Documentation

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

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

Definition at line 247 of file time_solver.h.

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

bool libMesh::ReferenceCounter::_enable_print_counter = true
staticprotectedinherited

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

Definition at line 143 of file reference_counter.h.

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

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

An implicit linear solver to use for adjoint problems.

Definition at line 252 of file time_solver.h.

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

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 137 of file reference_counter.h.

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

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

Definition at line 132 of file reference_counter.h.

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

UniquePtr<NumericVector<Number> > libMesh::SecondOrderUnsteadySolver::_old_local_solution_accel
protected

Serial vector of previous time step acceleration $ \ddot{u}_n $

Definition at line 112 of file second_order_unsteady_solver.h.

Referenced by libMesh::NewmarkSolver::advance_timestep(), init_data(), old_solution_accel(), and reinit().

UniquePtr<NumericVector<Number> > libMesh::SecondOrderUnsteadySolver::_old_local_solution_rate
protected

Serial vector of previous time step velocity $ \dot{u}_n $

Definition at line 107 of file second_order_unsteady_solver.h.

Referenced by libMesh::NewmarkSolver::advance_timestep(), init_data(), old_solution_rate(), and reinit().

sys_type& libMesh::TimeSolver::_system
protectedinherited

A reference to the system we are solving.

Definition at line 257 of file time_solver.h.

Referenced by libMesh::EulerSolver::_general_residual(), libMesh::Euler2Solver::_general_residual(), libMesh::SteadySolver::_general_residual(), libMesh::NewmarkSolver::_general_residual(), libMesh::UnsteadySolver::adjoint_advance_timestep(), libMesh::NewmarkSolver::advance_timestep(), libMesh::AdaptiveTimeSolver::advance_timestep(), libMesh::UnsteadySolver::advance_timestep(), libMesh::NewmarkSolver::compute_initial_accel(), libMesh::FirstOrderUnsteadySolver::compute_second_order_eqns(), libMesh::UnsteadySolver::du(), libMesh::EulerSolver::element_residual(), libMesh::Euler2Solver::element_residual(), libMesh::EigenTimeSolver::element_residual(), init(), libMesh::UnsteadySolver::init(), libMesh::TimeSolver::init(), libMesh::EigenTimeSolver::init(), 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(), old_solution_accel(), old_solution_rate(), libMesh::NewmarkSolver::project_initial_accel(), project_initial_rate(), 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().

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

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

bool libMesh::UnsteadySolver::first_solve
protectedinherited
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
This has no effect for SteadySolvers.
You must set at least one of the DiffSolver flags "continue_after_max_iterations" or "continue_after_backtrack_failure" to allow the TimeSolver to retry the solve.

Definition at line 220 of file time_solver.h.

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

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

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

Definition at line 264 of file time_solver.h.

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


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