GRINS-0.8.0
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
GRINS::UnsteadySolver Class Reference

#include <unsteady_solver.h>

Inheritance diagram for GRINS::UnsteadySolver:
Inheritance graph
[legend]
Collaboration diagram for GRINS::UnsteadySolver:
Collaboration graph
[legend]

Public Member Functions

 UnsteadySolver (const GetPot &input)
 
virtual ~UnsteadySolver ()
 
virtual void solve (SolverContext &context)
 
- Public Member Functions inherited from GRINS::Solver
 Solver (const GetPot &input)
 
virtual ~Solver ()
 
virtual void initialize (const GetPot &input, SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system)
 
virtual void adjoint_qoi_parameter_sensitivity (SolverContext &, const libMesh::QoISet &, const libMesh::ParameterVector &, libMesh::SensitivityData &) const
 
virtual void forward_qoi_parameter_sensitivity (SolverContext &, const libMesh::QoISet &, const libMesh::ParameterVector &, libMesh::SensitivityData &) const
 
void steady_adjoint_solve (SolverContext &context)
 Do steady version of adjoint solve. More...
 
void print_scalar_vars (SolverContext &context)
 
void print_qoi (SolverContext &context)
 

Protected Member Functions

virtual void init_time_solver (GRINS::MultiphysicsSystem *system)
 
template<typename T >
void set_theta (libMesh::UnsteadySolver *time_solver)
 
void update_dirichlet_bcs (SolverContext &context)
 Updates Dirichlet boundary conditions. More...
 
void init_second_order_in_time_solvers (SolverContext &context)
 
- Protected Member Functions inherited from GRINS::Solver
void set_solver_options (libMesh::DiffSolver &solver)
 

Protected Attributes

std::string _time_solver_name
 
unsigned int _n_timesteps
 
unsigned int _backtrack_deltat
 
double _theta
 
double _deltat
 
AdaptiveTimeSteppingOptions _adapt_time_step_options
 
bool _is_second_order_in_time
 Track whether is this a second order (in time) solver or not. More...
 
- Protected Attributes inherited from GRINS::Solver
unsigned int _max_nonlinear_iterations
 
double _relative_step_tolerance
 
double _absolute_step_tolerance
 
double _relative_residual_tolerance
 
double _absolute_residual_tolerance
 
double _initial_linear_tolerance
 
double _minimum_linear_tolerance
 
unsigned int _max_linear_iterations
 
bool _continue_after_backtrack_failure
 
bool _continue_after_max_iterations
 
bool _require_residual_reduction
 
bool _solver_quiet
 
bool _solver_verbose
 

Detailed Description

Definition at line 39 of file unsteady_solver.h.

Constructor & Destructor Documentation

GRINS::UnsteadySolver::UnsteadySolver ( const GetPot &  input)

Definition at line 54 of file unsteady_solver.C.

55  : Solver(input),
63  {}
std::string _time_solver_name
static unsigned int parse_n_timesteps(const GetPot &input)
unsigned int _n_timesteps
Solver(const GetPot &input)
Definition: solver.C:48
static unsigned int parse_backtrack_deltat(const GetPot &input)
Parse option to retry failed time steps with smaller .
static std::string parse_time_stepper_name(const GetPot &input)
AdaptiveTimeSteppingOptions _adapt_time_step_options
bool _is_second_order_in_time
Track whether is this a second order (in time) solver or not.
unsigned int _backtrack_deltat
static double parse_theta(const GetPot &input)
Parse value of for theta method time stepping.
static double parse_deltat(const GetPot &input)
virtual GRINS::UnsteadySolver::~UnsteadySolver ( )
inlinevirtual

Definition at line 44 of file unsteady_solver.h.

44 {};

Member Function Documentation

void GRINS::UnsteadySolver::init_second_order_in_time_solvers ( SolverContext context)
protected

Definition at line 227 of file unsteady_solver.C.

References GRINS::SolverContext::have_restart, and GRINS::SolverContext::system.

Referenced by solve().

228  {
229  // Right now, only Newmark is available so we cast directly to that
230  libMesh::TimeSolver& base_time_solver = context.system->get_time_solver();
231 
232  libMesh::NewmarkSolver& time_solver = libMesh::cast_ref<libMesh::NewmarkSolver&>(base_time_solver);
233 
234  // If there's a restart, the acceleration should already be there
235  if( context.have_restart )
236  time_solver.set_initial_accel_avail(true);
237 
238  // Otherwise, we need to compute it
239  else
240  {
241  libMesh::out << "==========================================================" << std::endl
242  << " Computing Initital Acceleration" << std::endl
243  << "==========================================================" << std::endl;
244 
245  time_solver.compute_initial_accel();
246  }
247  }
void GRINS::UnsteadySolver::init_time_solver ( GRINS::MultiphysicsSystem system)
protectedvirtual

Implements GRINS::Solver.

Definition at line 65 of file unsteady_solver.C.

References _adapt_time_step_options, _backtrack_deltat, _is_second_order_in_time, _time_solver_name, GRINS::AdaptiveTimeSteppingOptions::component_norm(), GRINS::AdaptiveTimeSteppingOptions::is_time_adaptive(), GRINS::SolverNames::libmesh_euler2_solver(), GRINS::SolverNames::libmesh_euler_solver(), GRINS::SolverNames::libmesh_newmark_solver(), GRINS::AdaptiveTimeSteppingOptions::max_growth(), GRINS::AdaptiveTimeSteppingOptions::target_tolerance(), and GRINS::AdaptiveTimeSteppingOptions::upper_tolerance().

66  {
67  libMesh::UnsteadySolver* time_solver = NULL;
68 
70  {
71  time_solver = new libMesh::EulerSolver( *(system) );
72 
73  this->set_theta<libMesh::EulerSolver>(time_solver);
74  }
76  {
77  time_solver = new libMesh::Euler2Solver( *(system) );
78 
79  this->set_theta<libMesh::Euler2Solver>(time_solver);
80  }
82  {
83  time_solver = new libMesh::NewmarkSolver( *(system) );
85  }
86  else
87  libmesh_error_msg("ERROR: Unsupported time stepper "+_time_solver_name);
88 
90  {
91  libMesh::TwostepTimeSolver *outer_solver =
92  new libMesh::TwostepTimeSolver(*system);
93 
94  outer_solver->target_tolerance = _adapt_time_step_options.target_tolerance();
95  outer_solver->upper_tolerance = _adapt_time_step_options.upper_tolerance();
96  outer_solver->max_growth = _adapt_time_step_options.max_growth();
97  outer_solver->component_norm = _adapt_time_step_options.component_norm();
98  outer_solver->quiet = false;
99 
100  outer_solver->core_time_solver =
101  libMesh::UniquePtr<libMesh::UnsteadySolver>(time_solver);
102  system->time_solver = libMesh::UniquePtr<libMesh::TimeSolver>(outer_solver);
103  }
104  else
105  {
106  system->time_solver = libMesh::UniquePtr<libMesh::TimeSolver>(time_solver);
107  }
108 
109  time_solver->reduce_deltat_on_diffsolver_failure = this->_backtrack_deltat;
110  }
std::string _time_solver_name
static const std::string libmesh_euler2_solver()
Definition: solver_names.h:49
AdaptiveTimeSteppingOptions _adapt_time_step_options
bool _is_second_order_in_time
Track whether is this a second order (in time) solver or not.
unsigned int _backtrack_deltat
static const std::string libmesh_euler_solver()
Definition: solver_names.h:46
static const std::string libmesh_newmark_solver()
Definition: solver_names.h:52
template<typename T >
void GRINS::UnsteadySolver::set_theta ( libMesh::UnsteadySolver *  time_solver)
inlineprotected

Definition at line 80 of file unsteady_solver.h.

References _theta.

81  {
82  T* deriv_solver = libMesh::cast_ptr<T*>(time_solver);
83  deriv_solver->theta = this->_theta;
84  }
void GRINS::UnsteadySolver::solve ( SolverContext context)
virtual

Implements GRINS::Solver.

Reimplemented in GRINS::UnsteadyMeshAdaptiveSolver.

Definition at line 112 of file unsteady_solver.C.

References _deltat, _is_second_order_in_time, _n_timesteps, GRINS::SolverContext::equation_system, init_second_order_in_time_solvers(), GRINS::SolverContext::output_residual, GRINS::SolverContext::output_vis, GRINS::SolverContext::postprocessing, GRINS::SolverContext::print_perflog, GRINS::Solver::print_scalar_vars(), GRINS::SolverContext::print_scalars, GRINS::SolverContext::system, GRINS::SolverContext::timesteps_per_perflog, GRINS::SolverContext::timesteps_per_vis, update_dirichlet_bcs(), and GRINS::SolverContext::vis.

113  {
114  libmesh_assert( context.system );
115 
116  context.system->deltat = this->_deltat;
117 
118  libMesh::Real sim_time;
119 
120  if( context.output_vis )
121  {
122  context.postprocessing->update_quantities( *(context.equation_system) );
123  context.vis->output( context.equation_system );
124  }
125 
126  // We may need to initialize acceleration for second order solvers
128  this->init_second_order_in_time_solvers(context);
129 
130  std::time_t first_wall_time = std::time(NULL);
131 
132  // Now we begin the timestep loop to compute the time-accurate
133  // solution of the equations.
134  for (unsigned int t_step=0; t_step < this->_n_timesteps; t_step++)
135  {
136  std::time_t latest_wall_time = std::time(NULL);
137 
138  std::cout << "==========================================================" << std::endl
139  << " Beginning time step " << t_step <<
140  ", t = " << context.system->time <<
141  ", dt = " << context.system->deltat <<
142  ", runtime = " << (latest_wall_time - first_wall_time) <<
143  std::endl
144  << "==========================================================" << std::endl;
145 
146  // If we have any solution-dependent Dirichlet boundaries, we
147  // need to update them with the current solution.
148  this->update_dirichlet_bcs(context);
149 
150  context.system->solve();
151 
152  sim_time = context.system->time;
153 
154  if( context.output_vis && !((t_step+1)%context.timesteps_per_vis) )
155  {
156  context.postprocessing->update_quantities( *(context.equation_system) );
157  context.vis->output( context.equation_system, t_step, sim_time );
158  }
159 
160  if( context.output_residual && !((t_step+1)%context.timesteps_per_vis) )
161  context.vis->output_residual( context.equation_system, context.system,
162  t_step, sim_time );
163 
164  if ( context.print_perflog && context.timesteps_per_perflog
165  && !((t_step+1)%context.timesteps_per_perflog) )
166  libMesh::perflog.print_log();
167 
168  if ( context.print_scalars )
169  this->print_scalar_vars(context);
170 
171  // Advance to the next timestep
172  context.system->time_solver->advance_timestep();
173  }
174 
175  std::time_t final_wall_time = std::time(NULL);
176  std::cout << "==========================================================" << std::endl
177  << " Ending time stepping, t = " << context.system->time <<
178  ", runtime = " << (final_wall_time - first_wall_time) <<
179  std::endl
180  << "==========================================================" << std::endl;
181 
182 
183  return;
184  }
unsigned int _n_timesteps
bool _is_second_order_in_time
Track whether is this a second order (in time) solver or not.
void print_scalar_vars(SolverContext &context)
Definition: solver.C:133
void update_dirichlet_bcs(SolverContext &context)
Updates Dirichlet boundary conditions.
void init_second_order_in_time_solvers(SolverContext &context)
void GRINS::UnsteadySolver::update_dirichlet_bcs ( SolverContext context)
protected

Updates Dirichlet boundary conditions.

If the Dirichlet boundary condition is nonlinear or time-dependent, we need to update the constraints with the new solution.

Definition at line 186 of file unsteady_solver.C.

References GRINS::SolverContext::system.

Referenced by solve(), and GRINS::UnsteadyMeshAdaptiveSolver::solve().

187  {
188  // FIXME: This needs to be much more efficient and intuitive.
189  bool have_nonlinear_dirichlet_bc = false;
190  bool have_time_dependence = false;
191  {
192  const libMesh::DirichletBoundaries &db =
193  *context.system->get_dof_map().get_dirichlet_boundaries();
194 
195  for (libMesh::DirichletBoundaries::const_iterator
196  it = db.begin(); it != db.end(); ++it)
197  {
198  const libMesh::DirichletBoundary* bdy = *it;
199 
200  // If we have a FEMFunctionBase, we assume nonlinearity
201  if (bdy->f_fem.get())
202  have_nonlinear_dirichlet_bc = true;
203 
204  // Check for time-dependence of FunctionBase
205  if( bdy->f.get() )
206  if( bdy->f->is_time_dependent() )
207  have_time_dependence = true;
208 
209  if( have_nonlinear_dirichlet_bc || have_time_dependence )
210  break;
211 
212  } // End loop over DirichletBoundaries
213  }
214 
215 
216  // Nonlinear Dirichlet constraints change as the solution does
217  // and time-dependent constraints have to be updated
218  if (have_nonlinear_dirichlet_bc || have_time_dependence )
219  {
220  context.system->reinit_constraints();
221  context.system->get_dof_map().enforce_constraints_exactly(*context.system);
222  context.system->get_dof_map().enforce_constraints_exactly(*context.system,
223  dynamic_cast<libMesh::UnsteadySolver*>(context.system->time_solver.get())->old_local_nonlinear_solution.get());
224  }
225  }

Member Data Documentation

AdaptiveTimeSteppingOptions GRINS::UnsteadySolver::_adapt_time_step_options
protected

Definition at line 70 of file unsteady_solver.h.

Referenced by init_time_solver().

unsigned int GRINS::UnsteadySolver::_backtrack_deltat
protected

Definition at line 65 of file unsteady_solver.h.

Referenced by init_time_solver().

double GRINS::UnsteadySolver::_deltat
protected

Definition at line 67 of file unsteady_solver.h.

Referenced by solve(), and GRINS::UnsteadyMeshAdaptiveSolver::solve().

bool GRINS::UnsteadySolver::_is_second_order_in_time
protected

Track whether is this a second order (in time) solver or not.

If it is, we need to potentially initialize the acceleration

Definition at line 74 of file unsteady_solver.h.

Referenced by init_time_solver(), and solve().

unsigned int GRINS::UnsteadySolver::_n_timesteps
protected

Definition at line 64 of file unsteady_solver.h.

Referenced by solve(), and GRINS::UnsteadyMeshAdaptiveSolver::solve().

double GRINS::UnsteadySolver::_theta
protected

Definition at line 66 of file unsteady_solver.h.

Referenced by set_theta().

std::string GRINS::UnsteadySolver::_time_solver_name
protected

Definition at line 62 of file unsteady_solver.h.

Referenced by init_time_solver().


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

Generated on Tue Dec 19 2017 12:47:32 for GRINS-0.8.0 by  doxygen 1.8.9.1