GRINS-0.7.0
List of all members | Public Member Functions
GRINS::UnsteadyMeshAdaptiveSolver Class Reference

#include <unsteady_mesh_adaptive_solver.h>

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

Public Member Functions

 UnsteadyMeshAdaptiveSolver (const GetPot &input)
 
virtual ~UnsteadyMeshAdaptiveSolver ()
 
virtual void solve (SolverContext &context)
 
- Public Member Functions inherited from GRINS::UnsteadySolver
 UnsteadySolver (const GetPot &input)
 
virtual ~UnsteadySolver ()
 
- 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, std::ostream &output)
 
- Public Member Functions inherited from GRINS::MeshAdaptiveSolverBase
 MeshAdaptiveSolverBase (const GetPot &input)
 
virtual ~MeshAdaptiveSolverBase ()
 

Additional Inherited Members

- Protected Types inherited from GRINS::MeshAdaptiveSolverBase
enum  RefinementFlaggingType {
  INVALID = 0, ERROR_TOLERANCE, N_ELEM_TARGET, ERROR_FRACTION,
  ELEM_FRACTION, MEAN_STD_DEV
}
 
- Protected Member Functions inherited from GRINS::UnsteadySolver
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 Member Functions inherited from GRINS::MeshAdaptiveSolverBase
void build_mesh_refinement (libMesh::MeshBase &mesh)
 
void set_refinement_type (const GetPot &input, const MeshAdaptivityOptions &mesh_adaptivity_options, RefinementFlaggingType &refinement_type)
 
bool check_for_convergence (SolverContext &context, const libMesh::ErrorVector &error) const
 
void flag_elements_for_refinement (const libMesh::ErrorVector &error)
 
void estimate_error_for_amr (SolverContext &context, libMesh::ErrorVector &error)
 
void perform_amr (SolverContext &context, const libMesh::ErrorVector &error)
 
- Protected Attributes inherited from GRINS::UnsteadySolver
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
 
- Protected Attributes inherited from GRINS::MeshAdaptiveSolverBase
ErrorEstimatorOptions _error_estimator_options
 
MeshAdaptivityOptions _mesh_adaptivity_options
 
RefinementFlaggingType _refinement_type
 
libMesh::UniquePtr< libMesh::MeshRefinement > _mesh_refinement
 

Detailed Description

Definition at line 38 of file unsteady_mesh_adaptive_solver.h.

Constructor & Destructor Documentation

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

Definition at line 40 of file unsteady_mesh_adaptive_solver.C.

41  : UnsteadySolver(input),
42  MeshAdaptiveSolverBase( input )
43  {}
UnsteadySolver(const GetPot &input)
virtual GRINS::UnsteadyMeshAdaptiveSolver::~UnsteadyMeshAdaptiveSolver ( )
inlinevirtual

Definition at line 45 of file unsteady_mesh_adaptive_solver.h.

45 {};

Member Function Documentation

void GRINS::UnsteadyMeshAdaptiveSolver::solve ( SolverContext context)
virtual

Reimplemented from GRINS::UnsteadySolver.

Definition at line 45 of file unsteady_mesh_adaptive_solver.C.

References GRINS::UnsteadySolver::_deltat, GRINS::MeshAdaptiveSolverBase::_mesh_adaptivity_options, GRINS::UnsteadySolver::_n_timesteps, GRINS::MeshAdaptiveSolverBase::build_mesh_refinement(), GRINS::MeshAdaptiveSolverBase::check_for_convergence(), GRINS::SolverContext::equation_system, GRINS::MeshAdaptiveSolverBase::estimate_error_for_amr(), GRINS::MeshAdaptivityOptions::max_refinement_steps(), GRINS::SolverContext::output_residual, GRINS::SolverContext::output_vis, GRINS::MeshAdaptiveSolverBase::perform_amr(), GRINS::SolverContext::postprocessing, GRINS::SolverContext::print_perflog, GRINS::SolverContext::print_qoi, GRINS::Solver::print_qoi(), GRINS::Solver::print_scalar_vars(), GRINS::SolverContext::print_scalars, GRINS::SolverContext::system, GRINS::SolverContext::timesteps_per_perflog, GRINS::SolverContext::timesteps_per_vis, GRINS::UnsteadySolver::update_dirichlet_bcs(), and GRINS::SolverContext::vis.

46  {
47  context.system->deltat = this->_deltat;
48 
49  libMesh::Real sim_time;
50 
51  if( context.output_vis )
52  {
53  context.postprocessing->update_quantities( *(context.equation_system) );
54  context.vis->output( context.equation_system );
55  }
56 
57  // Setup MeshRefinement
58  libMesh::MeshBase& mesh = context.equation_system->get_mesh();
59  this->build_mesh_refinement( mesh );
60 
61  std::time_t first_wall_time = std::time(NULL);
62 
63  // Now we begin the timestep loop to compute the time-accurate
64  // solution of the equations.
65  for (unsigned int t_step=0; t_step < this->_n_timesteps; t_step++)
66  {
67  std::time_t latest_wall_time = std::time(NULL);
68 
69  std::cout << "==========================================================" << std::endl
70  << " Beginning time step " << t_step <<
71  ", t = " << context.system->time <<
72  ", dt = " << context.system->deltat <<
73  ", runtime = " << (latest_wall_time - first_wall_time) <<
74  std::endl
75  << "==========================================================" << std::endl;
76 
77  // If we have any solution-dependent Dirichlet boundaries, we
78  // need to update them with the current solution.
79  this->update_dirichlet_bcs(context);
80 
81  for ( unsigned int r_step = 0; r_step < _mesh_adaptivity_options.max_refinement_steps(); r_step++ )
82  {
83  std::cout << "==========================================================" << std::endl
84  << "Adaptive Refinement Step " << r_step << std::endl
85  << "==========================================================" << std::endl;
86 
87  // GRVY timers contained in here (if enabled)
88  context.system->solve();
89 
90  sim_time = context.system->time;
91 
92  if( context.output_vis && !((t_step+1)%context.timesteps_per_vis) )
93  {
94  context.postprocessing->update_quantities( *(context.equation_system) );
95  context.vis->output( context.equation_system, t_step, sim_time );
96  }
97 
98  if( context.output_residual && !((t_step+1)%context.timesteps_per_vis) )
99  context.vis->output_residual( context.equation_system, context.system,
100  t_step, sim_time );
101 
102  if ( context.print_perflog && context.timesteps_per_perflog
103  && !((t_step+1)%context.timesteps_per_perflog) )
104  libMesh::perflog.print_log();
105 
106  if ( context.print_scalars )
107  this->print_scalar_vars(context);
108 
109 
110  // Now we construct the data structures for the mesh refinement process
111  libMesh::ErrorVector error;
112  this->estimate_error_for_amr(context,error);
113 
114  // Check for convergence of error
115  bool converged = this->check_for_convergence( context, error );
116 
117  if( converged )
118  {
119  // Break out of adaptive loop
120  std::cout << "==========================================================" << std::endl
121  << "Convergence detected!" << std::endl
122  << "==========================================================" << std::endl;
123  break;
124  }
125  else
126  {
127  // Only bother refining if we're not on the last step.
129  this->perform_amr(context,error);
130  }
131 
132  } // End mesh adaptive loop
133 
134  // Advance to the next timestep
135  context.system->time_solver->advance_timestep();
136 
137  } // End time step loop
138 
139  std::time_t final_wall_time = std::time(NULL);
140  std::cout << "==========================================================" << std::endl
141  << " Ending time stepping, t = " << context.system->time <<
142  ", runtime = " << (final_wall_time - first_wall_time) <<
143  std::endl
144  << "==========================================================" << std::endl;
145 
146  // Print out the QoI, but only do it if the user asks for it
147  if( context.print_qoi )
148  this->print_qoi(context,std::cout);
149  }
void build_mesh_refinement(libMesh::MeshBase &mesh)
void perform_amr(SolverContext &context, const libMesh::ErrorVector &error)
void estimate_error_for_amr(SolverContext &context, libMesh::ErrorVector &error)
bool check_for_convergence(SolverContext &context, const libMesh::ErrorVector &error) const
void print_scalar_vars(SolverContext &context)
Definition: grins_solver.C:133
void print_qoi(SolverContext &context, std::ostream &output)
Definition: grins_solver.C:155
unsigned int max_refinement_steps() const
void update_dirichlet_bcs(SolverContext &context)
Updates Dirichlet boundary conditions.
MeshAdaptivityOptions _mesh_adaptivity_options

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

Generated on Thu Jun 2 2016 21:52:32 for GRINS-0.7.0 by  doxygen 1.8.10