EulerSolver Class Reference
#include <euler_solver.h>

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 | solve () |
| virtual void | advance_timestep () |
| Number | old_nonlinear_solution (const unsigned int global_dof_number) const |
| virtual Real | du (const SystemNorm &norm) const |
| virtual void | reinit () |
| virtual void | adjoint_recede_timestep () |
| virtual void | before_timestep () |
| const sys_type & | system () const |
| virtual AutoPtr< DiffSolver > & | diff_solver () |
| virtual AutoPtr< LinearSolver < Number > > & | linear_solver () |
Static Public Member Functions | |
| static std::string | get_info () |
| static void | print_info () |
| static unsigned int | n_objects () |
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 | |
| sys_type & | system () |
| void | increment_constructor_count (const std::string &name) |
| void | increment_destructor_count (const std::string &name) |
Protected Attributes | |
| bool | first_solve |
| AutoPtr< DiffSolver > | _diff_solver |
| AutoPtr< LinearSolver< Number > > | _linear_solver |
| sys_type & | _system |
Static Protected Attributes | |
| static Counts | _counts |
| static Threads::atomic < unsigned int > | _n_objects |
| static Threads::spin_mutex | _mutex |
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.
Definition at line 44 of file euler_solver.h.
Member Typedef Documentation
typedef std::map<std::string, std::pair<unsigned int, unsigned int> > ReferenceCounter::Counts [protected, inherited] |
Data structure to log the information. The log is identified by the class name.
Definition at line 105 of file reference_counter.h.
| typedef UnsteadySolver EulerSolver::Parent |
The parent class
Definition at line 50 of file euler_solver.h.
typedef DifferentiableSystem TimeSolver::sys_type [inherited] |
The type of system
Reimplemented in EigenTimeSolver, and SteadySolver.
Definition at line 65 of file time_solver.h.
Constructor & Destructor Documentation
| EulerSolver::EulerSolver | ( | sys_type & | s | ) |
Constructor. Requires a reference to the system to be solved.
Definition at line 27 of file euler_solver.C.
00028 : UnsteadySolver(s), theta(1.) 00029 { 00030 }
| EulerSolver::~EulerSolver | ( | ) | [virtual] |
Member Function Documentation
| void TimeSolver::adjoint_recede_timestep | ( | ) | [virtual, inherited] |
This method advances the adjoint solution to the previous timestep, after an adjoint_solve() has been performed. This will probably be done after every UnsteadySolver::adjoint_solve().
Definition at line 84 of file time_solver.C.
| void UnsteadySolver::advance_timestep | ( | ) | [virtual, inherited] |
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 TimeSolver.
Reimplemented in AdaptiveTimeSolver.
Definition at line 124 of file unsteady_solver.C.
References TimeSolver::_system, DifferentiableSystem::deltat, UnsteadySolver::first_solve, System::get_vector(), UnsteadySolver::old_nonlinear_solution(), System::solution, and DifferentiableSystem::time.
Referenced by UnsteadySolver::solve().
00125 { 00126 NumericVector<Number> &old_nonlinear_solution = 00127 _system.get_vector("_old_nonlinear_solution"); 00128 NumericVector<Number> &nonlinear_solution = 00129 *(_system.solution); 00130 00131 old_nonlinear_solution = nonlinear_solution; 00132 00133 if (!first_solve) 00134 _system.time += _system.deltat; 00135 }
| virtual void TimeSolver::before_timestep | ( | ) | [inline, virtual, inherited] |
This method is for subclasses or users to override to do arbitrary processing between timesteps
Definition at line 137 of file time_solver.h.
| virtual AutoPtr<DiffSolver>& TimeSolver::diff_solver | ( | ) | [inline, virtual, inherited] |
An implicit linear or nonlinear solver to use at each timestep.
Reimplemented in AdaptiveTimeSolver.
Definition at line 147 of file time_solver.h.
References TimeSolver::_diff_solver.
00147 { return _diff_solver; }
| Real UnsteadySolver::du | ( | const SystemNorm & | norm | ) | const [virtual, inherited] |
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 TimeSolver.
Definition at line 150 of file unsteady_solver.C.
References TimeSolver::_system, System::calculate_norm(), System::get_vector(), and System::solution.
00151 { 00152 00153 AutoPtr<NumericVector<Number> > solution_copy = 00154 _system.solution->clone(); 00155 00156 solution_copy->add(-1., _system.get_vector("_old_nonlinear_solution")); 00157 00158 solution_copy->close(); 00159 00160 return _system.calculate_norm(*solution_copy, norm); 00161 }
| bool 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 TimeSolver.
Definition at line 50 of file euler_solver.C.
References TimeSolver::_system, DenseVector< T >::add(), DifferentiableSystem::deltat, DiffContext::dof_indices, DiffContext::elem_fixed_solution, DiffContext::elem_jacobian, DiffContext::elem_reinit(), DiffContext::elem_residual, DiffContext::elem_solution, DiffContext::elem_solution_derivative, DifferentiableSystem::element_constraint(), DifferentiableSystem::element_time_derivative(), DifferentiableSystem::eulerian_residual(), DiffContext::fixed_solution_derivative, DifferentiableSystem::mass_residual(), UnsteadySolver::old_nonlinear_solution(), theta, and DifferentiableSystem::use_fixed_solution.
00052 { 00053 unsigned int n_dofs = context.elem_solution.size(); 00054 00055 // Local nonlinear solution at old timestep 00056 DenseVector<Number> old_elem_solution(n_dofs); 00057 for (unsigned int i=0; i != n_dofs; ++i) 00058 old_elem_solution(i) = 00059 old_nonlinear_solution(context.dof_indices[i]); 00060 00061 // Local nonlinear solution at time t_theta 00062 DenseVector<Number> theta_solution(context.elem_solution); 00063 theta_solution *= theta; 00064 theta_solution.add(1. - theta, old_elem_solution); 00065 00066 // Technically the elem_solution_derivative is either theta 00067 // or -1.0 in this implementation, but we scale the former part 00068 // ourselves 00069 context.elem_solution_derivative = 1.0; 00070 00071 // Technically the fixed_solution_derivative is always theta, 00072 // but we're scaling a whole jacobian by theta after these first 00073 // evaluations 00074 context.fixed_solution_derivative = 1.0; 00075 00076 // If a fixed solution is requested, we'll use theta_solution 00077 if (_system.use_fixed_solution) 00078 context.elem_fixed_solution = theta_solution; 00079 00080 // Move theta_->elem_, elem_->theta_ 00081 context.elem_solution.swap(theta_solution); 00082 00083 // Move the mesh into place first if necessary 00084 context.elem_reinit(theta); 00085 00086 // We're going to compute just the change in elem_residual 00087 // (and possibly elem_jacobian), then add back the old values 00088 DenseVector<Number> old_elem_residual(context.elem_residual); 00089 DenseMatrix<Number> old_elem_jacobian; 00090 if (request_jacobian) 00091 { 00092 old_elem_jacobian = context.elem_jacobian; 00093 context.elem_jacobian.zero(); 00094 } 00095 context.elem_residual.zero(); 00096 00097 // Get the time derivative at t_theta 00098 bool jacobian_computed = 00099 _system.element_time_derivative(request_jacobian, context); 00100 00101 // For a moving mesh problem we may need the pseudoconvection term too 00102 jacobian_computed = 00103 _system.eulerian_residual(jacobian_computed, context) && jacobian_computed; 00104 00105 // Scale the time-dependent residual and jacobian correctly 00106 context.elem_residual *= _system.deltat; 00107 if (jacobian_computed) 00108 context.elem_jacobian *= (theta * _system.deltat); 00109 00110 // The fixed_solution_derivative is always theta, 00111 // and now we're done scaling jacobians 00112 context.fixed_solution_derivative = theta; 00113 00114 // We evaluate mass_residual with the change in solution 00115 // to get the mass matrix, reusing old_elem_solution to hold that 00116 // delta_solution. We're solving dt*F(u) - du = 0, so our 00117 // delta_solution is old_solution - new_solution. 00118 // We're still keeping elem_solution in theta_solution for now 00119 old_elem_solution -= theta_solution; 00120 00121 // Move old_->elem_, theta_->old_ 00122 context.elem_solution.swap(old_elem_solution); 00123 00124 // We do a trick here to avoid using a non-1 00125 // elem_solution_derivative: 00126 context.elem_jacobian *= -1.0; 00127 jacobian_computed = _system.mass_residual(jacobian_computed, context) && 00128 jacobian_computed; 00129 context.elem_jacobian *= -1.0; 00130 00131 // Move elem_->elem_, old_->theta_ 00132 context.elem_solution.swap(theta_solution); 00133 00134 // Restore the elem position if necessary 00135 context.elem_reinit(1.); 00136 00137 // Add the constraint term 00138 jacobian_computed = _system.element_constraint(jacobian_computed, context) && 00139 jacobian_computed; 00140 00141 // Add back the old residual and jacobian 00142 context.elem_residual += old_elem_residual; 00143 if (request_jacobian) 00144 { 00145 if (jacobian_computed) 00146 context.elem_jacobian += old_elem_jacobian; 00147 else 00148 context.elem_jacobian.swap(old_elem_jacobian); 00149 } 00150 00151 return jacobian_computed; 00152 }
| Real EulerSolver::error_order | ( | ) | const [virtual] |
Error convergence order: 2 for Crank-Nicolson, 1 otherwise
Implements UnsteadySolver.
Definition at line 40 of file euler_solver.C.
References theta.
00041 { 00042 if (theta == 0.5) 00043 return 2.; 00044 return 1.; 00045 }
| std::string ReferenceCounter::get_info | ( | ) | [static, inherited] |
Gets a string containing the reference information.
Definition at line 45 of file reference_counter.C.
References ReferenceCounter::_counts, and QuadratureRules::name().
Referenced by ReferenceCounter::print_info().
00046 { 00047 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00048 00049 std::ostringstream out; 00050 00051 out << '\n' 00052 << " ---------------------------------------------------------------------------- \n" 00053 << "| Reference count information |\n" 00054 << " ---------------------------------------------------------------------------- \n"; 00055 00056 for (Counts::iterator it = _counts.begin(); 00057 it != _counts.end(); ++it) 00058 { 00059 const std::string name(it->first); 00060 const unsigned int creations = it->second.first; 00061 const unsigned int destructions = it->second.second; 00062 00063 out << "| " << name << " reference count information:\n" 00064 << "| Creations: " << creations << '\n' 00065 << "| Destructions: " << destructions << '\n'; 00066 } 00067 00068 out << " ---------------------------------------------------------------------------- \n"; 00069 00070 return out.str(); 00071 00072 #else 00073 00074 return ""; 00075 00076 #endif 00077 }
| void ReferenceCounter::increment_constructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.
Definition at line 149 of file reference_counter.h.
References ReferenceCounter::_counts, and Threads::spin_mtx.
Referenced by ReferenceCountedObject< SparseMatrix< T > >::ReferenceCountedObject().
00150 { 00151 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00152 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00153 00154 p.first++; 00155 }
| void ReferenceCounter::increment_destructor_count | ( | const std::string & | name | ) | [inline, protected, inherited] |
Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.
Definition at line 167 of file reference_counter.h.
References ReferenceCounter::_counts, and Threads::spin_mtx.
Referenced by ReferenceCountedObject< SparseMatrix< T > >::~ReferenceCountedObject().
00168 { 00169 Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx); 00170 std::pair<unsigned int, unsigned int>& p = _counts[name]; 00171 00172 p.second++; 00173 }
| void UnsteadySolver::init | ( | ) | [virtual, inherited] |
The initialization function. This method is used to initialize internal data structures before a simulation begins.
Reimplemented from TimeSolver.
Reimplemented in AdaptiveTimeSolver.
Definition at line 44 of file unsteady_solver.C.
References TimeSolver::_system, and System::add_vector().
00045 { 00046 TimeSolver::init(); 00047 00048 _system.add_vector("_old_nonlinear_solution"); 00049 }
| virtual AutoPtr<LinearSolver<Number> >& TimeSolver::linear_solver | ( | ) | [inline, virtual, inherited] |
An implicit linear solver to use for adjoint and sensitivity problems.
Definition at line 152 of file time_solver.h.
References TimeSolver::_linear_solver.
00152 { return _linear_solver; }
| static unsigned int ReferenceCounter::n_objects | ( | ) | [inline, static, inherited] |
Prints the number of outstanding (created, but not yet destroyed) objects.
Definition at line 76 of file reference_counter.h.
References ReferenceCounter::_n_objects.
00077 { return _n_objects; }
| Number UnsteadySolver::old_nonlinear_solution | ( | const unsigned int | global_dof_number | ) | const [inherited] |
- Returns:
- the old nonlinear solution for the specified global DOF.
Definition at line 139 of file unsteady_solver.C.
References TimeSolver::_system, System::get_dof_map(), DofMap::n_dofs(), and UnsteadySolver::old_local_nonlinear_solution.
Referenced by UnsteadySolver::advance_timestep(), AdaptiveTimeSolver::advance_timestep(), element_residual(), Euler2Solver::element_residual(), FEMSystem::eulerian_residual(), side_residual(), and Euler2Solver::side_residual().
00141 { 00142 libmesh_assert (global_dof_number < _system.get_dof_map().n_dofs()); 00143 libmesh_assert (global_dof_number < old_local_nonlinear_solution->size()); 00144 00145 return (*old_local_nonlinear_solution)(global_dof_number); 00146 }
| void ReferenceCounter::print_info | ( | ) | [static, inherited] |
Prints the reference information to std::cout.
Definition at line 83 of file reference_counter.C.
References ReferenceCounter::get_info().
00084 { 00085 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG) 00086 00087 std::cout << ReferenceCounter::get_info(); 00088 00089 #endif 00090 }
| void TimeSolver::reinit | ( | ) | [virtual, inherited] |
The reinitialization function. This method is used after changes in the mesh
Reimplemented in AdaptiveTimeSolver, and EigenTimeSolver.
Definition at line 46 of file time_solver.C.
References TimeSolver::_diff_solver, and TimeSolver::_linear_solver.
00047 { 00048 _diff_solver->reinit(); 00049 _linear_solver->clear(); 00050 _linear_solver->init(); 00051 }
| bool 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 TimeSolver.
Definition at line 156 of file euler_solver.C.
References TimeSolver::_system, DenseVector< T >::add(), DifferentiableSystem::deltat, DiffContext::dof_indices, DiffContext::elem_fixed_solution, DiffContext::elem_jacobian, DiffContext::elem_residual, DiffContext::elem_side_reinit(), DiffContext::elem_solution, DiffContext::elem_solution_derivative, DiffContext::fixed_solution_derivative, UnsteadySolver::old_nonlinear_solution(), DifferentiableSystem::side_constraint(), DifferentiableSystem::side_mass_residual(), DifferentiableSystem::side_time_derivative(), theta, and DifferentiableSystem::use_fixed_solution.
00158 { 00159 unsigned int n_dofs = context.elem_solution.size(); 00160 00161 // Local nonlinear solution at old timestep 00162 DenseVector<Number> old_elem_solution(n_dofs); 00163 for (unsigned int i=0; i != n_dofs; ++i) 00164 old_elem_solution(i) = 00165 old_nonlinear_solution(context.dof_indices[i]); 00166 00167 // Local nonlinear solution at time t_theta 00168 DenseVector<Number> theta_solution(context.elem_solution); 00169 theta_solution *= theta; 00170 theta_solution.add(1. - theta, old_elem_solution); 00171 00172 // Technically the elem_solution_derivative is either theta 00173 // or 1.0 in this implementation, but we scale the former part 00174 // ourselves 00175 context.elem_solution_derivative = 1.0; 00176 00177 // Technically the fixed_solution_derivative is always theta, 00178 // but we're scaling a whole jacobian by theta after these first 00179 // evaluations 00180 context.fixed_solution_derivative = 1.0; 00181 00182 // If a fixed solution is requested, we'll use theta_solution 00183 if (_system.use_fixed_solution) 00184 context.elem_fixed_solution = theta_solution; 00185 00186 // Move theta_->elem_, elem_->theta_ 00187 context.elem_solution.swap(theta_solution); 00188 00189 // Move the mesh into place first if necessary 00190 context.elem_side_reinit(theta); 00191 00192 // We're going to compute just the change in elem_residual 00193 // (and possibly elem_jacobian), then add back the old values 00194 DenseVector<Number> old_elem_residual(context.elem_residual); 00195 DenseMatrix<Number> old_elem_jacobian; 00196 if (request_jacobian) 00197 { 00198 old_elem_jacobian = context.elem_jacobian; 00199 context.elem_jacobian.zero(); 00200 } 00201 context.elem_residual.zero(); 00202 00203 // Get the time derivative at t_theta 00204 bool jacobian_computed = 00205 _system.side_time_derivative(request_jacobian, context); 00206 00207 // Scale the time-dependent residual and jacobian correctly 00208 context.elem_residual *= _system.deltat; 00209 if (jacobian_computed) 00210 context.elem_jacobian *= (theta * _system.deltat); 00211 00212 // The fixed_solution_derivative is always theta, 00213 // and now we're done scaling jacobians 00214 context.fixed_solution_derivative = theta; 00215 00216 // We evaluate side_mass_residual with the change in solution 00217 // to get the mass matrix, reusing old_elem_solution to hold that 00218 // delta_solution. We're solving dt*F(u) - du = 0, so our 00219 // delta_solution is old_solution - new_solution. 00220 // We're still keeping elem_solution in theta_solution for now 00221 old_elem_solution -= theta_solution; 00222 00223 // Move old_->elem_, theta_->old_ 00224 context.elem_solution.swap(old_elem_solution); 00225 00226 // We do a trick here to avoid using a non-1 00227 // elem_solution_derivative: 00228 context.elem_jacobian *= -1.0; 00229 jacobian_computed = _system.side_mass_residual(jacobian_computed, context) && 00230 jacobian_computed; 00231 context.elem_jacobian *= -1.0; 00232 00233 // Move elem_->elem_, old_->theta_ 00234 context.elem_solution.swap(theta_solution); 00235 00236 // Restore the elem position if necessary 00237 context.elem_side_reinit(1.); 00238 00239 // Add the constraint term 00240 jacobian_computed = _system.side_constraint(jacobian_computed, context) && 00241 jacobian_computed; 00242 00243 // Add back the old residual and jacobian 00244 context.elem_residual += old_elem_residual; 00245 if (request_jacobian) 00246 { 00247 if (jacobian_computed) 00248 context.elem_jacobian += old_elem_jacobian; 00249 else 00250 context.elem_jacobian.swap(old_elem_jacobian); 00251 } 00252 00253 return jacobian_computed; 00254 }
| void UnsteadySolver::solve | ( | ) | [virtual, inherited] |
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 TimeSolver.
Reimplemented in AdaptiveTimeSolver, and TwostepTimeSolver.
Definition at line 53 of file unsteady_solver.C.
References TimeSolver::_diff_solver, TimeSolver::_system, UnsteadySolver::advance_timestep(), DifferentiableSystem::deltat, DiffSolver::DIVERGED_BACKTRACKING_FAILURE, DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS, UnsteadySolver::first_solve, System::get_dof_map(), DofMap::get_send_list(), System::get_vector(), libMeshEnums::GHOSTED, System::n_dofs(), System::n_local_dofs(), UnsteadySolver::old_local_nonlinear_solution, TimeSolver::quiet, TimeSolver::reduce_deltat_on_diffsolver_failure, and libMeshEnums::SERIAL.
00054 { 00055 if (first_solve) 00056 { 00057 advance_timestep(); 00058 first_solve = false; 00059 } 00060 00061 #ifdef LIBMESH_ENABLE_GHOSTED 00062 old_local_nonlinear_solution->init (_system.n_dofs(), _system.n_local_dofs(), 00063 _system.get_dof_map().get_send_list(), false, 00064 GHOSTED); 00065 #else 00066 old_local_nonlinear_solution->init (_system.n_dofs(), false, SERIAL); 00067 #endif 00068 00069 _system.get_vector("_old_nonlinear_solution").localize 00070 (*old_local_nonlinear_solution, 00071 _system.get_dof_map().get_send_list()); 00072 00073 unsigned int solve_result = _diff_solver->solve(); 00074 00075 // If we requested the UnsteadySolver to attempt reducing dt after a 00076 // failed DiffSolver solve, check the results of the solve now. 00077 if (reduce_deltat_on_diffsolver_failure) 00078 { 00079 bool backtracking_failed = 00080 solve_result & DiffSolver::DIVERGED_BACKTRACKING_FAILURE; 00081 00082 bool max_iterations = 00083 solve_result & DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS; 00084 00085 if (backtracking_failed || max_iterations) 00086 { 00087 // Cut timestep in half 00088 for (unsigned int nr=0; nr<reduce_deltat_on_diffsolver_failure; ++nr) 00089 { 00090 _system.deltat *= 0.5; 00091 std::cout << "Newton backtracking failed. Trying with smaller timestep, dt=" 00092 << _system.deltat << std::endl; 00093 00094 solve_result = _diff_solver->solve(); 00095 00096 // Check solve results with reduced timestep 00097 bool backtracking_failed = 00098 solve_result & DiffSolver::DIVERGED_BACKTRACKING_FAILURE; 00099 00100 bool max_iterations = 00101 solve_result & DiffSolver::DIVERGED_MAX_NONLINEAR_ITERATIONS; 00102 00103 if (!backtracking_failed && !max_iterations) 00104 { 00105 if (!quiet) 00106 std::cout << "Reduced dt solve succeeded." << std::endl; 00107 return; 00108 } 00109 } 00110 00111 // If we made it here, we still couldn't converge the solve after 00112 // reducing deltat 00113 std::cout << "DiffSolver::solve() did not succeed after " 00114 << reduce_deltat_on_diffsolver_failure 00115 << " attempts." << std::endl; 00116 libmesh_convergence_failure(); 00117 00118 } // end if (backtracking_failed || max_iterations) 00119 } // end if (reduce_deltat_on_diffsolver_failure) 00120 }
| sys_type& TimeSolver::system | ( | ) | [inline, protected, inherited] |
- Returns:
- a writeable reference to the system we are solving.
Definition at line 197 of file time_solver.h.
References TimeSolver::_system.
00197 { return _system; }
| const sys_type& TimeSolver::system | ( | ) | const [inline, inherited] |
- Returns:
- a constant reference to the system we are solving.
Definition at line 142 of file time_solver.h.
References TimeSolver::_system.
00142 { return _system; }
Member Data Documentation
ReferenceCounter::Counts ReferenceCounter::_counts [static, protected, inherited] |
Actually holds the data.
Definition at line 110 of file reference_counter.h.
Referenced by ReferenceCounter::get_info(), ReferenceCounter::increment_constructor_count(), and ReferenceCounter::increment_destructor_count().
AutoPtr<DiffSolver> TimeSolver::_diff_solver [protected, inherited] |
An implicit linear or nonlinear solver to use at each timestep.
Definition at line 187 of file time_solver.h.
Referenced by TimeSolver::diff_solver(), TimeSolver::init(), TimeSolver::reinit(), UnsteadySolver::solve(), and TimeSolver::solve().
AutoPtr<LinearSolver<Number> > TimeSolver::_linear_solver [protected, inherited] |
An implicit linear solver to use for adjoint problems.
Definition at line 192 of file time_solver.h.
Referenced by TimeSolver::init(), TimeSolver::linear_solver(), and TimeSolver::reinit().
Threads::spin_mutex ReferenceCounter::_mutex [static, protected, inherited] |
Mutual exclusion object to enable thread-safe reference counting.
Definition at line 123 of file reference_counter.h.
Threads::atomic< unsigned int > ReferenceCounter::_n_objects [static, protected, inherited] |
The number of objects. Print the reference count information when the number returns to 0.
Definition at line 118 of file reference_counter.h.
Referenced by ReferenceCounter::n_objects(), ReferenceCounter::ReferenceCounter(), and ReferenceCounter::~ReferenceCounter().
sys_type& TimeSolver::_system [protected, inherited] |
A reference to the system we are solving.
Definition at line 202 of file time_solver.h.
Referenced by UnsteadySolver::advance_timestep(), AdaptiveTimeSolver::advance_timestep(), UnsteadySolver::du(), SteadySolver::element_residual(), element_residual(), Euler2Solver::element_residual(), EigenTimeSolver::element_residual(), UnsteadySolver::init(), TimeSolver::init(), EigenTimeSolver::init(), UnsteadySolver::old_nonlinear_solution(), SteadySolver::side_residual(), side_residual(), Euler2Solver::side_residual(), EigenTimeSolver::side_residual(), UnsteadySolver::solve(), TwostepTimeSolver::solve(), EigenTimeSolver::solve(), and TimeSolver::system().
bool UnsteadySolver::first_solve [protected, inherited] |
A bool that will be true the first time solve() is called, and false thereafter
Reimplemented from TimeSolver.
Definition at line 124 of file unsteady_solver.h.
Referenced by UnsteadySolver::advance_timestep(), AdaptiveTimeSolver::advance_timestep(), UnsteadySolver::solve(), and TwostepTimeSolver::solve().
AutoPtr<NumericVector<Number> > UnsteadySolver::old_local_nonlinear_solution [inherited] |
Serial vector of _system.get_vector("_old_nonlinear_solution")
Reimplemented from TimeSolver.
Definition at line 105 of file unsteady_solver.h.
Referenced by AdaptiveTimeSolver::AdaptiveTimeSolver(), AdaptiveTimeSolver::init(), UnsteadySolver::old_nonlinear_solution(), UnsteadySolver::solve(), and AdaptiveTimeSolver::~AdaptiveTimeSolver().
bool TimeSolver::quiet [inherited] |
Print extra debugging information if quiet == false.
Definition at line 157 of file time_solver.h.
Referenced by UnsteadySolver::solve(), TwostepTimeSolver::solve(), and EigenTimeSolver::solve().
unsigned int 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 180 of file time_solver.h.
Referenced by UnsteadySolver::solve(), and TwostepTimeSolver::solve().
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 91 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: