libMesh::EulerSolver Class Reference

#include <euler_solver.h>

Inheritance diagram for libMesh::EulerSolver:

Public Types

typedef UnsteadySolver Parent
 
typedef DifferentiableSystem sys_type
 

Public Member Functions

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

Protected Types

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

bool first_solve
 
bool first_adjoint_step
 
AutoPtr< DiffSolver_diff_solver
 
AutoPtr< LinearSolver< Number > > _linear_solver
 
sys_type_system
 
AutoPtr< 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 2006

Definition at line 45 of file euler_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 113 of file reference_counter.h.

The parent class

Definition at line 51 of file euler_solver.h.

The type of system

Definition at line 66 of file time_solver.h.

Constructor & Destructor Documentation

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.

28  : UnsteadySolver(s), theta(1.)
29 {
30 }
libMesh::EulerSolver::~EulerSolver ( )
virtual

Destructor.

Definition at line 34 of file euler_solver.C.

35 {
36 }

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.

Definition at line 176 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.

177 {
178  // On the first call of this function, we dont save the adjoint solution or
179  // decrement the time, we just call the retrieve function below
180  if(!first_adjoint_step)
181  {
182  // Call the store function to store the last adjoint before decrementing the time
183  solution_history->store();
184  // Decrement the system time
186  }
187  else
188  {
189  first_adjoint_step = false;
190  }
191 
192  // Retrieve the primal solution vectors at this time using the
193  // solution_history object
194  solution_history->retrieve();
195 
196  // Dont forget to localize the old_nonlinear_solution !
197  _system.get_vector("_old_nonlinear_solution").localize
200 }
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.

Definition at line 150 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::UnsteadySolver::solve().

151 {
152  if (!first_solve)
153  {
154  // Store the solution, does nothing by default
155  // User has to attach appropriate solution_history object for this to
156  // actually store anything anywhere
157  solution_history->store();
158 
160  }
161 
162  NumericVector<Number> &old_nonlinear_soln =
163  _system.get_vector("_old_nonlinear_solution");
164  NumericVector<Number> &nonlinear_solution =
165  *(_system.solution);
166 
167  old_nonlinear_soln = nonlinear_solution;
168 
169  old_nonlinear_soln.localize
172 }
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 152 of file time_solver.h.

152 {}
virtual AutoPtr<DiffSolver>& libMesh::TimeSolver::diff_solver ( )
inlinevirtualinherited

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

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 167 of file time_solver.h.

References libMesh::TimeSolver::_diff_solver.

167 { return _diff_solver; }
void libMesh::ReferenceCounter::disable_print_counter_info ( )
staticinherited

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

107 {
108  _enable_print_counter = false;
109  return;
110 }
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 225 of file unsteady_solver.C.

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

226 {
227 
228  AutoPtr<NumericVector<Number> > solution_copy =
229  _system.solution->clone();
230 
231  solution_copy->add(-1., _system.get_vector("_old_nonlinear_solution"));
232 
233  solution_copy->close();
234 
235  return _system.calculate_norm(*solution_copy, norm);
236 }
bool libMesh::EulerSolver::element_residual ( bool  request_jacobian,
DiffContext context 
)
virtual

This method uses the DifferentiableSystem's 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::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, libMesh::DiffContext::elem_reinit(), libMesh::DiffContext::elem_solution_derivative, libMesh::DifferentiablePhysics::element_constraint(), libMesh::DifferentiablePhysics::element_time_derivative(), libMesh::DifferentiablePhysics::eulerian_residual(), 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_residual(), libMesh::DiffContext::get_elem_solution(), libMesh::DifferentiablePhysics::mass_residual(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::DenseVector< T >::size(), libMesh::DenseVector< T >::swap(), libMesh::DenseMatrix< T >::swap(), theta, libMesh::System::use_fixed_solution, libMesh::DenseMatrix< T >::zero(), and libMesh::DenseVector< T >::zero().

52 {
53  unsigned int n_dofs = context.get_elem_solution().size();
54 
55  // Local nonlinear solution at old timestep
56  DenseVector<Number> old_elem_solution(n_dofs);
57  for (unsigned int i=0; i != n_dofs; ++i)
58  old_elem_solution(i) =
59  old_nonlinear_solution(context.get_dof_indices()[i]);
60 
61  // Local nonlinear solution at time t_theta
62  DenseVector<Number> theta_solution(context.get_elem_solution());
63  theta_solution *= theta;
64  theta_solution.add(1. - theta, old_elem_solution);
65 
66  // Technically the elem_solution_derivative is either theta
67  // or -1.0 in this implementation, but we scale the former part
68  // ourselves
69  context.elem_solution_derivative = 1.0;
70 
71 // Technically the fixed_solution_derivative is always theta,
72 // but we're scaling a whole jacobian by theta after these first
73 // evaluations
74  context.fixed_solution_derivative = 1.0;
75 
76  // If a fixed solution is requested, we'll use theta_solution
78  context.get_elem_fixed_solution() = theta_solution;
79 
80  // Move theta_->elem_, elem_->theta_
81  context.get_elem_solution().swap(theta_solution);
82 
83  // Move the mesh into place first if necessary
84  context.elem_reinit(theta);
85 
86  // We're going to compute just the change in elem_residual
87  // (and possibly elem_jacobian), then add back the old values
88  DenseVector<Number> old_elem_residual(context.get_elem_residual());
89  DenseMatrix<Number> old_elem_jacobian;
90  if (request_jacobian)
91  {
92  old_elem_jacobian = context.get_elem_jacobian();
93  context.get_elem_jacobian().zero();
94  }
95  context.get_elem_residual().zero();
96 
97  // Get the time derivative at t_theta
98  bool jacobian_computed =
99  _system.element_time_derivative(request_jacobian, context);
100 
101  // For a moving mesh problem we may need the pseudoconvection term too
102  jacobian_computed =
103  _system.eulerian_residual(jacobian_computed, context) && jacobian_computed;
104 
105  // Scale the time-dependent residual and jacobian correctly
106  context.get_elem_residual() *= _system.deltat;
107  if (jacobian_computed)
108  context.get_elem_jacobian() *= (theta * _system.deltat);
109 
110  // The fixed_solution_derivative is always theta,
111  // and now we're done scaling jacobians
112  context.fixed_solution_derivative = theta;
113 
114  // We evaluate mass_residual with the change in solution
115  // to get the mass matrix, reusing old_elem_solution to hold that
116  // delta_solution. We're solving dt*F(u) - du = 0, so our
117  // delta_solution is old_solution - new_solution.
118  // We're still keeping elem_solution in theta_solution for now
119  old_elem_solution -= theta_solution;
120 
121  // Move old_->elem_, theta_->old_
122  context.get_elem_solution().swap(old_elem_solution);
123 
124  // We do a trick here to avoid using a non-1
125  // elem_solution_derivative:
126  context.get_elem_jacobian() *= -1.0;
127  context.fixed_solution_derivative *= -1.0;
128  jacobian_computed = _system.mass_residual(jacobian_computed, context) &&
129  jacobian_computed;
130  context.get_elem_jacobian() *= -1.0;
131  context.fixed_solution_derivative *= -1.0;
132 
133  // Move elem_->elem_, old_->theta_
134  context.get_elem_solution().swap(theta_solution);
135 
136  // Restore the elem position if necessary
137  context.elem_reinit(1.);
138 
139  // Add the constraint term
140  jacobian_computed = _system.element_constraint(jacobian_computed, context) &&
141  jacobian_computed;
142 
143  // Add back the old residual and jacobian
144  context.get_elem_residual() += old_elem_residual;
145  if (request_jacobian)
146  {
147  if (jacobian_computed)
148  context.get_elem_jacobian() += old_elem_jacobian;
149  else
150  context.get_elem_jacobian().swap(old_elem_jacobian);
151  }
152 
153  return jacobian_computed;
154 }
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 }
Real libMesh::EulerSolver::error_order ( ) const
virtual

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

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

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

164 {
165  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
166  std::pair<unsigned int, unsigned int>& p = _counts[name];
167 
168  p.first++;
169 }
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 176 of file reference_counter.h.

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

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

177 {
178  Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
179  std::pair<unsigned int, unsigned int>& p = _counts[name];
180 
181  p.second++;
182 }
void libMesh::UnsteadySolver::init ( )
virtualinherited

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

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 46 of file unsteady_solver.C.

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

47 {
49 
50  _system.add_vector("_old_nonlinear_solution");
51 }
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.

Definition at line 55 of file unsteady_solver.C.

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

56 {
57 #ifdef LIBMESH_ENABLE_GHOSTED
60  GHOSTED);
61 #else
63 #endif
64 }
bool libMesh::TimeSolver::is_adjoint ( ) const
inlineinherited

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

Definition at line 217 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

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

This is not a steady-state solver.

Implements libMesh::TimeSolver.

Definition at line 149 of file unsteady_solver.h.

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

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

Definition at line 172 of file time_solver.h.

References libMesh::TimeSolver::_linear_solver.

172 { return _linear_solver; }
static unsigned int libMesh::ReferenceCounter::n_objects ( )
inlinestaticinherited

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

Definition at line 79 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

80  { return _n_objects; }
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 214 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 element_residual(), libMesh::Euler2Solver::element_residual(), side_residual(), and libMesh::Euler2Solver::side_residual().

216 {
217  libmesh_assert_less (global_dof_number, _system.get_dof_map().n_dofs());
218  libmesh_assert_less (global_dof_number, old_local_nonlinear_solution->size());
219 
220  return (*old_local_nonlinear_solution)(global_dof_number);
221 }
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().

89 {
91 }
void libMesh::UnsteadySolver::reinit ( )
virtualinherited

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

Reimplemented from libMesh::TimeSolver.

Reimplemented in libMesh::AdaptiveTimeSolver.

Definition at line 68 of file unsteady_solver.C.

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

69 {
71 
72 #ifdef LIBMESH_ENABLE_GHOSTED
75  GHOSTED);
76 #else
78 #endif
79 
80  // localize the old solution
81  NumericVector<Number> &old_nonlinear_soln =
82  _system.get_vector("_old_nonlinear_solution");
83 
84  old_nonlinear_soln.localize
87 }
void libMesh::UnsteadySolver::retrieve_timestep ( )
virtualinherited

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

Reimplemented from libMesh::TimeSolver.

Definition at line 202 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.

203  {
204  // Retrieve all the stored vectors at the current time
205  solution_history->retrieve();
206 
207  // Dont forget to localize the old_nonlinear_solution !
208  _system.get_vector("_old_nonlinear_solution").localize
211  }
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 224 of file time_solver.h.

References libMesh::TimeSolver::_is_adjoint.

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

225  { _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 91 of file time_solver.C.

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

92  {
93  solution_history = _solution_history.clone();
94  }
bool libMesh::EulerSolver::side_residual ( bool  request_jacobian,
DiffContext context 
)
virtual

This method uses the DifferentiableSystem's 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 158 of file euler_solver.C.

References libMesh::TimeSolver::_system, libMesh::DifferentiableSystem::deltat, libMesh::DiffContext::elem_side_reinit(), libMesh::DiffContext::elem_solution_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_residual(), libMesh::DiffContext::get_elem_solution(), libMesh::UnsteadySolver::old_nonlinear_solution(), libMesh::DifferentiablePhysics::side_constraint(), libMesh::DifferentiablePhysics::side_mass_residual(), libMesh::DifferentiablePhysics::side_time_derivative(), libMesh::DenseVector< T >::size(), libMesh::DenseVector< T >::swap(), libMesh::DenseMatrix< T >::swap(), theta, libMesh::System::use_fixed_solution, libMesh::DenseMatrix< T >::zero(), and libMesh::DenseVector< T >::zero().

160 {
161  unsigned int n_dofs = context.get_elem_solution().size();
162 
163  // Local nonlinear solution at old timestep
164  DenseVector<Number> old_elem_solution(n_dofs);
165  for (unsigned int i=0; i != n_dofs; ++i)
166  old_elem_solution(i) =
167  old_nonlinear_solution(context.get_dof_indices()[i]);
168 
169  // Local nonlinear solution at time t_theta
170  DenseVector<Number> theta_solution(context.get_elem_solution());
171  theta_solution *= theta;
172  theta_solution.add(1. - theta, old_elem_solution);
173 
174  // Technically the elem_solution_derivative is either theta
175  // or 1.0 in this implementation, but we scale the former part
176  // ourselves
177  context.elem_solution_derivative = 1.0;
178 
179 // Technically the fixed_solution_derivative is always theta,
180 // but we're scaling a whole jacobian by theta after these first
181 // evaluations
182  context.fixed_solution_derivative = 1.0;
183 
184  // If a fixed solution is requested, we'll use theta_solution
186  context.get_elem_fixed_solution() = theta_solution;
187 
188  // Move theta_->elem_, elem_->theta_
189  context.get_elem_solution().swap(theta_solution);
190 
191  // Move the mesh into place first if necessary
192  context.elem_side_reinit(theta);
193 
194  // We're going to compute just the change in elem_residual
195  // (and possibly elem_jacobian), then add back the old values
196  DenseVector<Number> old_elem_residual(context.get_elem_residual());
197  DenseMatrix<Number> old_elem_jacobian;
198  if (request_jacobian)
199  {
200  old_elem_jacobian = context.get_elem_jacobian();
201  context.get_elem_jacobian().zero();
202  }
203  context.get_elem_residual().zero();
204 
205  // Get the time derivative at t_theta
206  bool jacobian_computed =
207  _system.side_time_derivative(request_jacobian, context);
208 
209  // Scale the time-dependent residual and jacobian correctly
210  context.get_elem_residual() *= _system.deltat;
211  if (jacobian_computed)
212  context.get_elem_jacobian() *= (theta * _system.deltat);
213 
214  // The fixed_solution_derivative is always theta,
215  // and now we're done scaling jacobians
216  context.fixed_solution_derivative = theta;
217 
218  // We evaluate side_mass_residual with the change in solution
219  // to get the mass matrix, reusing old_elem_solution to hold that
220  // delta_solution. We're solving dt*F(u) - du = 0, so our
221  // delta_solution is old_solution - new_solution.
222  // We're still keeping elem_solution in theta_solution for now
223  old_elem_solution -= theta_solution;
224 
225  // Move old_->elem_, theta_->old_
226  context.get_elem_solution().swap(old_elem_solution);
227 
228  // We do a trick here to avoid using a non-1
229  // elem_solution_derivative:
230  context.get_elem_jacobian() *= -1.0;
231  jacobian_computed = _system.side_mass_residual(jacobian_computed, context) &&
232  jacobian_computed;
233  context.get_elem_jacobian() *= -1.0;
234 
235  // Move elem_->elem_, old_->theta_
236  context.get_elem_solution().swap(theta_solution);
237 
238  // Restore the elem position if necessary
239  context.elem_side_reinit(1.);
240 
241  // Add the constraint term
242  jacobian_computed = _system.side_constraint(jacobian_computed, context) &&
243  jacobian_computed;
244 
245  // Add back the old residual and jacobian
246  context.get_elem_residual() += old_elem_residual;
247  if (request_jacobian)
248  {
249  if (jacobian_computed)
250  context.get_elem_jacobian() += old_elem_jacobian;
251  else
252  context.get_elem_jacobian().swap(old_elem_jacobian);
253  }
254 
255  return jacobian_computed;
256 }
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::AdaptiveTimeSolver, and libMesh::TwostepTimeSolver.

Definition at line 91 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.

92 {
93  if (first_solve)
94  {
96  first_solve = false;
97  }
98 
99  unsigned int solve_result = _diff_solver->solve();
100 
101  // If we requested the UnsteadySolver to attempt reducing dt after a
102  // failed DiffSolver solve, check the results of the solve now.
104  {
105  bool backtracking_failed =
107 
108  bool max_iterations =
110 
111  if (backtracking_failed || max_iterations)
112  {
113  // Cut timestep in half
114  for (unsigned int nr=0; nr<reduce_deltat_on_diffsolver_failure; ++nr)
115  {
116  _system.deltat *= 0.5;
117  libMesh::out << "Newton backtracking failed. Trying with smaller timestep, dt="
118  << _system.deltat << std::endl;
119 
120  solve_result = _diff_solver->solve();
121 
122  // Check solve results with reduced timestep
123  bool backtracking_still_failed =
125 
126  bool backtracking_max_iterations =
128 
129  if (!backtracking_still_failed && !backtracking_max_iterations)
130  {
131  if (!quiet)
132  libMesh::out << "Reduced dt solve succeeded." << std::endl;
133  return;
134  }
135  }
136 
137  // If we made it here, we still couldn't converge the solve after
138  // reducing deltat
139  libMesh::out << "DiffSolver::solve() did not succeed after "
140  << reduce_deltat_on_diffsolver_failure
141  << " attempts." << std::endl;
142  libmesh_convergence_failure();
143 
144  } // end if (backtracking_failed || max_iterations)
145  } // end if (reduce_deltat_on_diffsolver_failure)
146 }
const sys_type& libMesh::TimeSolver::system ( ) const
inlineinherited
Returns
a constant reference to the system we are solving.

Definition at line 157 of file time_solver.h.

References libMesh::TimeSolver::_system.

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

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

Definition at line 162 of file time_solver.h.

References libMesh::TimeSolver::_system.

162 { return _system; }

Member Data Documentation

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

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

Definition at line 232 of file time_solver.h.

Referenced by libMesh::TimeSolver::diff_solver(), libMesh::TimeSolver::init(), libMesh::TimeSolver::reinit(), libMesh::UnsteadySolver::solve(), and libMesh::TimeSolver::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 137 of file reference_counter.h.

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

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

An implicit linear solver to use for adjoint problems.

Definition at line 237 of file time_solver.h.

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

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

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 131 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 126 of file reference_counter.h.

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

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

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

bool libMesh::TimeSolver::quiet
inherited

Print extra debugging information if quiet == false.

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

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

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

An AutoPtr 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 260 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::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 93 of file euler_solver.h.

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


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

Site Created By: libMesh Developers
Last modified: February 07 2014 16:58:00 UTC

Hosted By:
SourceForge.net Logo