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

#include <grins_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)
 
virtual void initialize (const GetPot &input, std::tr1::shared_ptr< 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...
 

Protected Member Functions

virtual void init_time_solver (GRINS::MultiphysicsSystem *system)
 
void set_solver_options (libMesh::DiffSolver &solver)
 

Protected Attributes

unsigned int _n_timesteps
 
unsigned int _backtrack_deltat
 
double _theta
 
double _deltat
 
double _target_tolerance
 
double _upper_tolerance
 
double _max_growth
 
libMesh::SystemNorm _component_norm
 
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 _solver_quiet
 
bool _solver_verbose
 
std::map< std::string, GRINS::NBCContainer_neumann_bc_funcs
 

Detailed Description

Definition at line 37 of file grins_unsteady_solver.h.

Constructor & Destructor Documentation

GRINS::UnsteadySolver::UnsteadySolver ( const GetPot &  input)
Todo:
Is this the best default for delta t?

Definition at line 47 of file grins_unsteady_solver.C.

References _component_norm.

48  : Solver(input),
49  _n_timesteps( input("unsteady-solver/n_timesteps", 1 ) ),
50  _backtrack_deltat( input("unsteady-solver/backtrack_deltat", 0 ) ),
51  _theta( input("unsteady-solver/theta", 0.5 ) ),
53  _deltat( input("unsteady-solver/deltat", 0.0 ) ),
54  _target_tolerance( input("unsteady-solver/target_tolerance", 0.0 ) ),
55  _upper_tolerance( input("unsteady-solver/upper_tolerance", 0.0 ) ),
56  _max_growth( input("unsteady-solver/max_growth", 0.0 ) )
57  {
58  const unsigned int n_component_norm =
59  input.vector_variable_size("unsteady-solver/component_norm");
60  for (unsigned int i=0; i != n_component_norm; ++i)
61  {
62  const std::string current_norm = input("component_norm", std::string("L2"), i);
63  // TODO: replace this with string_to_enum with newer libMesh
64  if (current_norm == "GRINSEnums::L2")
65  _component_norm.set_type(i, libMesh::L2);
66  else if (current_norm == "GRINSEnums::H1")
67  _component_norm.set_type(i, libMesh::H1);
68  else
69  libmesh_not_implemented();
70  }
71 
72 
73  }
Solver(const GetPot &input)
Definition: grins_solver.C:44
libMesh::SystemNorm _component_norm
GRINS::UnsteadySolver::~UnsteadySolver ( )
virtual

Definition at line 75 of file grins_unsteady_solver.C.

76  {
77  return;
78  }

Member Function Documentation

virtual void GRINS::Solver::adjoint_qoi_parameter_sensitivity ( SolverContext ,
const libMesh::QoISet &  ,
const libMesh::ParameterVector &  ,
libMesh::SensitivityData &   
) const
inlinevirtualinherited

Reimplemented in GRINS::SteadyMeshAdaptiveSolver, and GRINS::SteadySolver.

Definition at line 71 of file grins_solver.h.

76  { libmesh_not_implemented(); }
virtual void GRINS::Solver::forward_qoi_parameter_sensitivity ( SolverContext ,
const libMesh::QoISet &  ,
const libMesh::ParameterVector &  ,
libMesh::SensitivityData &   
) const
inlinevirtualinherited

Reimplemented in GRINS::SteadyMeshAdaptiveSolver, and GRINS::SteadySolver.

Definition at line 79 of file grins_solver.h.

84  { libmesh_not_implemented(); }
void GRINS::UnsteadySolver::init_time_solver ( GRINS::MultiphysicsSystem system)
protectedvirtual

Implements GRINS::Solver.

Definition at line 80 of file grins_unsteady_solver.C.

References _backtrack_deltat, _max_growth, _target_tolerance, _theta, and _upper_tolerance.

81  {
82  libMesh::EulerSolver* time_solver = new libMesh::EulerSolver( *(system) );
83 
85  {
86  libMesh::TwostepTimeSolver *outer_solver =
87  new libMesh::TwostepTimeSolver(*system);
88 
89  outer_solver->target_tolerance = _target_tolerance;
90  outer_solver->upper_tolerance = _upper_tolerance;
91  outer_solver->max_growth = _max_growth;
92  outer_solver->quiet = false;
93 
94  outer_solver->core_time_solver =
95  libMesh::AutoPtr<libMesh::UnsteadySolver>(time_solver);
96  system->time_solver = libMesh::AutoPtr<libMesh::TimeSolver>(outer_solver);
97  }
98  else
99  {
100  system->time_solver = libMesh::AutoPtr<libMesh::TimeSolver>(time_solver);
101  }
102 
103  // Set theta parameter for time-stepping scheme
104  time_solver->theta = this->_theta;
105  time_solver->reduce_deltat_on_diffsolver_failure = this->_backtrack_deltat;
106 
107  return;
108  }
void GRINS::Solver::initialize ( const GetPot &  input,
std::tr1::shared_ptr< libMesh::EquationSystems >  equation_system,
GRINS::MultiphysicsSystem system 
)
virtualinherited

Reimplemented in GRINS::DisplacementContinuationSolver.

Definition at line 67 of file grins_solver.C.

References GRINS::Solver::init_time_solver(), and GRINS::Solver::set_solver_options().

Referenced by GRINS::DisplacementContinuationSolver::initialize().

70  {
71 
72  // Defined in subclasses depending on the solver used.
73  this->init_time_solver(system);
74 
75  // Initialize the system
76  equation_system->init();
77 
78  // Get diff solver to set options
79  libMesh::DiffSolver &solver = *(system->time_solver->diff_solver().get());
80 
81  // Set linear/nonlinear solver options
82  this->set_solver_options( solver );
83 
84  return;
85  }
virtual void init_time_solver(GRINS::MultiphysicsSystem *system)=0
void set_solver_options(libMesh::DiffSolver &solver)
Definition: grins_solver.C:87
void GRINS::Solver::set_solver_options ( libMesh::DiffSolver &  solver)
protectedinherited

Definition at line 87 of file grins_solver.C.

References GRINS::Solver::_absolute_residual_tolerance, GRINS::Solver::_absolute_step_tolerance, GRINS::Solver::_continue_after_backtrack_failure, GRINS::Solver::_continue_after_max_iterations, GRINS::Solver::_initial_linear_tolerance, GRINS::Solver::_max_linear_iterations, GRINS::Solver::_max_nonlinear_iterations, GRINS::Solver::_minimum_linear_tolerance, GRINS::Solver::_relative_residual_tolerance, GRINS::Solver::_relative_step_tolerance, GRINS::Solver::_solver_quiet, and GRINS::Solver::_solver_verbose.

Referenced by GRINS::Solver::initialize().

88  {
89  solver.quiet = this->_solver_quiet;
90  solver.verbose = this->_solver_verbose;
91  solver.max_nonlinear_iterations = this->_max_nonlinear_iterations;
92  solver.relative_step_tolerance = this->_relative_step_tolerance;
93  solver.absolute_step_tolerance = this->_absolute_step_tolerance;
94  solver.relative_residual_tolerance = this->_relative_residual_tolerance;
95  solver.absolute_residual_tolerance = this->_absolute_residual_tolerance;
96  solver.max_linear_iterations = this->_max_linear_iterations;
97  solver.continue_after_backtrack_failure = this->_continue_after_backtrack_failure;
98  solver.initial_linear_tolerance = this->_initial_linear_tolerance;
99  solver.minimum_linear_tolerance = this->_minimum_linear_tolerance;
100  solver.continue_after_max_iterations = this->_continue_after_max_iterations;
101 
102  return;
103  }
double _absolute_residual_tolerance
Definition: grins_solver.h:103
bool _solver_verbose
Definition: grins_solver.h:112
double _relative_residual_tolerance
Definition: grins_solver.h:101
unsigned int _max_linear_iterations
Definition: grins_solver.h:106
double _relative_step_tolerance
Definition: grins_solver.h:95
bool _continue_after_backtrack_failure
Definition: grins_solver.h:107
double _minimum_linear_tolerance
Definition: grins_solver.h:105
double _absolute_step_tolerance
Definition: grins_solver.h:96
bool _continue_after_max_iterations
Definition: grins_solver.h:108
unsigned int _max_nonlinear_iterations
Definition: grins_solver.h:94
double _initial_linear_tolerance
Definition: grins_solver.h:104
void GRINS::UnsteadySolver::solve ( SolverContext context)
virtual

Implements GRINS::Solver.

Definition at line 110 of file grins_unsteady_solver.C.

References _deltat, _n_timesteps, GRINS::SolverContext::equation_system, GRINS::SolverContext::output_residual, GRINS::SolverContext::output_vis, GRINS::SolverContext::postprocessing, GRINS::SolverContext::print_perflog, GRINS::SolverContext::print_scalars, GRINS::SolverContext::system, GRINS::SolverContext::timesteps_per_perflog, GRINS::SolverContext::timesteps_per_vis, and GRINS::SolverContext::vis.

111  {
112  libmesh_assert( context.system );
113 
114  context.system->deltat = this->_deltat;
115 
116  libMesh::Real sim_time;
117 
118  if( context.output_vis )
119  {
120  context.postprocessing->update_quantities( *(context.equation_system) );
121  context.vis->output( context.equation_system );
122  }
123 
124  std::time_t first_wall_time = std::time(NULL);
125 
126  // Now we begin the timestep loop to compute the time-accurate
127  // solution of the equations.
128  for (unsigned int t_step=0; t_step < this->_n_timesteps; t_step++)
129  {
130  std::time_t latest_wall_time = std::time(NULL);
131 
132  std::cout << "==========================================================" << std::endl
133  << " Beginning time step " << t_step <<
134  ", t = " << context.system->time <<
135  ", dt = " << context.system->deltat <<
136  ", runtime = " << (latest_wall_time - first_wall_time) <<
137  std::endl
138  << "==========================================================" << std::endl;
139 
140  // If we have any solution-dependent Dirichlet boundaries, we
141  // need to update them with the current solution.
142  //
143  // FIXME: This needs to be much more efficient and intuitive.
144  bool have_nonlinear_dirichlet_bc = false;
145  {
146  const libMesh::DirichletBoundaries &db =
147  *context.system->get_dof_map().get_dirichlet_boundaries();
148  for (libMesh::DirichletBoundaries::const_iterator
149  it = db.begin(); it != db.end(); ++it)
150  {
151  const libMesh::DirichletBoundary* bdy = *it;
152  if (bdy->f_fem.get())
153  {
154  have_nonlinear_dirichlet_bc = true;
155  break;
156  }
157  }
158  }
159  if (have_nonlinear_dirichlet_bc)
160  context.system->get_equation_systems().reinit();
161 
162  // GRVY timers contained in here (if enabled)
163  context.system->solve();
164 
165  sim_time = context.system->time;
166 
167  if( context.output_vis && !((t_step+1)%context.timesteps_per_vis) )
168  {
169  context.postprocessing->update_quantities( *(context.equation_system) );
170  context.vis->output( context.equation_system, t_step, sim_time );
171  }
172 
173  if( context.output_residual && !((t_step+1)%context.timesteps_per_vis) )
174  context.vis->output_residual( context.equation_system, context.system,
175  t_step, sim_time );
176 
177  if ( context.print_perflog && context.timesteps_per_perflog
178  && !((t_step+1)%context.timesteps_per_perflog) )
179  libMesh::perflog.print_log();
180 
181  if ( context.print_scalars )
182  for (unsigned int v=0; v != context.system->n_vars(); ++v)
183  if (context.system->variable(v).type().family ==
184  libMesh::SCALAR)
185  {
186  std::cout << context.system->variable_name(v) <<
187  " = {";
188  std::vector<libMesh::dof_id_type> scalar_indices;
189  context.system->get_dof_map().SCALAR_dof_indices
190  (scalar_indices, v);
191  if (scalar_indices.size())
192  std::cout <<
193  context.system->current_solution(scalar_indices[0]);
194  for (unsigned int i=1; i < scalar_indices.size();
195  ++i)
196  std::cout << ", " <<
197  context.system->current_solution(scalar_indices[i]);
198  std::cout << '}' << std::endl;
199  }
200 
201  // Advance to the next timestep
202  context.system->time_solver->advance_timestep();
203  }
204 
205  std::time_t final_wall_time = std::time(NULL);
206  std::cout << "==========================================================" << std::endl
207  << " Ending time stepping, t = " << context.system->time <<
208  ", runtime = " << (final_wall_time - first_wall_time) <<
209  std::endl
210  << "==========================================================" << std::endl;
211 
212 
213  return;
214  }
void GRINS::Solver::steady_adjoint_solve ( SolverContext context)
inherited

Do steady version of adjoint solve.

We put this here since we may want to reuse this in multiple different steady solves.

Definition at line 105 of file grins_solver.C.

References GRINS::SolverContext::system.

Referenced by GRINS::SteadySolver::solve(), and GRINS::SteadyMeshAdaptiveSolver::solve().

106  {
107  libMesh::out << "==========================================================" << std::endl
108  << "Solving adjoint problem." << std::endl
109  << "==========================================================" << std::endl;
110 
111  context.system->adjoint_solve();
112  context.system->set_adjoint_already_solved(true);
113  }

Member Data Documentation

double GRINS::Solver::_absolute_residual_tolerance
protectedinherited

Definition at line 103 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

double GRINS::Solver::_absolute_step_tolerance
protectedinherited

Definition at line 96 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

unsigned int GRINS::UnsteadySolver::_backtrack_deltat
protected

Definition at line 51 of file grins_unsteady_solver.h.

Referenced by init_time_solver().

libMesh::SystemNorm GRINS::UnsteadySolver::_component_norm
protected

Definition at line 59 of file grins_unsteady_solver.h.

Referenced by UnsteadySolver().

bool GRINS::Solver::_continue_after_backtrack_failure
protectedinherited

Definition at line 107 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

bool GRINS::Solver::_continue_after_max_iterations
protectedinherited

Definition at line 108 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

double GRINS::UnsteadySolver::_deltat
protected

Definition at line 53 of file grins_unsteady_solver.h.

Referenced by solve().

double GRINS::Solver::_initial_linear_tolerance
protectedinherited

Definition at line 104 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

double GRINS::UnsteadySolver::_max_growth
protected

Definition at line 58 of file grins_unsteady_solver.h.

Referenced by init_time_solver().

unsigned int GRINS::Solver::_max_linear_iterations
protectedinherited

Definition at line 106 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

unsigned int GRINS::Solver::_max_nonlinear_iterations
protectedinherited

Definition at line 94 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

double GRINS::Solver::_minimum_linear_tolerance
protectedinherited

Definition at line 105 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

unsigned int GRINS::UnsteadySolver::_n_timesteps
protected

Definition at line 50 of file grins_unsteady_solver.h.

Referenced by solve().

std::map< std::string, GRINS::NBCContainer > GRINS::Solver::_neumann_bc_funcs
protectedinherited

Definition at line 117 of file grins_solver.h.

double GRINS::Solver::_relative_residual_tolerance
protectedinherited

Definition at line 101 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

double GRINS::Solver::_relative_step_tolerance
protectedinherited

Definition at line 95 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

bool GRINS::Solver::_solver_quiet
protectedinherited

Definition at line 111 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

bool GRINS::Solver::_solver_verbose
protectedinherited

Definition at line 112 of file grins_solver.h.

Referenced by GRINS::Solver::set_solver_options().

double GRINS::UnsteadySolver::_target_tolerance
protected

Definition at line 56 of file grins_unsteady_solver.h.

Referenced by init_time_solver().

double GRINS::UnsteadySolver::_theta
protected

Definition at line 52 of file grins_unsteady_solver.h.

Referenced by init_time_solver().

double GRINS::UnsteadySolver::_upper_tolerance
protected

Definition at line 57 of file grins_unsteady_solver.h.

Referenced by init_time_solver().


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

Generated on Mon Jun 22 2015 21:32:24 for GRINS-0.6.0 by  doxygen 1.8.9.1