GRINS-0.8.0
unsteady_solver.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2017 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 // This class
27 #include "grins/unsteady_solver.h"
28 
29 // GRINS
30 #include "grins/grins_enums.h"
31 #include "grins/solver_context.h"
32 #include "grins/multiphysics_sys.h"
35 #include "grins/solver_names.h"
36 
37 // libMesh
38 #include "libmesh/dirichlet_boundaries.h"
39 #include "libmesh/dof_map.h"
40 #include "libmesh/getpot.h"
41 #include "libmesh/euler_solver.h"
42 #include "libmesh/euler2_solver.h"
43 #include "libmesh/function_base.h"
44 #include "libmesh/twostep_time_solver.h"
45 #include "libmesh/newmark_solver.h"
46 #include "libmesh/function_base.h"
47 
48 // C++
49 #include <ctime>
50 
51 namespace GRINS
52 {
53 
54  UnsteadySolver::UnsteadySolver( const GetPot& input )
55  : Solver(input),
56  _time_solver_name(TimeSteppingParsing::parse_time_stepper_name(input)),
57  _n_timesteps( TimeSteppingParsing::parse_n_timesteps(input) ),
58  _backtrack_deltat( TimeSteppingParsing::parse_backtrack_deltat(input) ),
59  _theta( TimeSteppingParsing::parse_theta(input) ),
60  _deltat( TimeSteppingParsing::parse_deltat(input) ),
61  _adapt_time_step_options(input),
62  _is_second_order_in_time(false)
63  {}
64 
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  }
111 
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  }
185 
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  }
226 
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  }
248 
249 } // namespace GRINS
std::string _time_solver_name
unsigned int _n_timesteps
SharedPtr< libMesh::EquationSystems > equation_system
SharedPtr< PostProcessedQuantities< libMesh::Real > > postprocessing
UnsteadySolver(const GetPot &input)
static const std::string libmesh_euler2_solver()
Definition: solver_names.h:49
GRINS namespace.
unsigned int timesteps_per_perflog
virtual void init_time_solver(GRINS::MultiphysicsSystem *system)
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 timesteps_per_vis
unsigned int _backtrack_deltat
GRINS::MultiphysicsSystem * system
SharedPtr< GRINS::Visualization > vis
void print_scalar_vars(SolverContext &context)
Definition: solver.C:133
Interface with libMesh for solving Multiphysics problems.
virtual void solve(SolverContext &context)
static const std::string libmesh_euler_solver()
Definition: solver_names.h:46
Simple class to hold objects passed to Solver::solve.
static const std::string libmesh_newmark_solver()
Definition: solver_names.h:52
void update_dirichlet_bcs(SolverContext &context)
Updates Dirichlet boundary conditions.
void init_second_order_in_time_solvers(SolverContext &context)

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