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

#include <steady_mesh_adaptive_solver.h>

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

Public Member Functions

 SteadyMeshAdaptiveSolver (const GetPot &input)
 
virtual ~SteadyMeshAdaptiveSolver ()
 
virtual void solve (SolverContext &context)
 
virtual void adjoint_qoi_parameter_sensitivity (SolverContext &context, const libMesh::QoISet &qoi_indices, const libMesh::ParameterVector &parameters_in, libMesh::SensitivityData &sensitivities) const
 
virtual void forward_qoi_parameter_sensitivity (SolverContext &context, const libMesh::QoISet &qoi_indices, const libMesh::ParameterVector &parameters_in, libMesh::SensitivityData &sensitivities) const
 
- 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)
 
void steady_adjoint_solve (SolverContext &context)
 Do steady version of adjoint solve. More...
 
void print_scalar_vars (SolverContext &context)
 
void print_qoi (SolverContext &context)
 
- Public Member Functions inherited from GRINS::MeshAdaptiveSolverBase
 MeshAdaptiveSolverBase (const GetPot &input)
 
virtual ~MeshAdaptiveSolverBase ()
 

Protected Member Functions

virtual void init_time_solver (MultiphysicsSystem *system)
 
void check_qoi_error_option_consistency (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)
 

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 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 steady_mesh_adaptive_solver.h.

Constructor & Destructor Documentation

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

Definition at line 44 of file steady_mesh_adaptive_solver.C.

45  : Solver(input),
46  MeshAdaptiveSolverBase( input )
47  {
48  return;
49  }
Solver(const GetPot &input)
Definition: solver.C:48
GRINS::SteadyMeshAdaptiveSolver::~SteadyMeshAdaptiveSolver ( )
virtual

Definition at line 51 of file steady_mesh_adaptive_solver.C.

52  {
53  return;
54  }

Member Function Documentation

void GRINS::SteadyMeshAdaptiveSolver::adjoint_qoi_parameter_sensitivity ( SolverContext context,
const libMesh::QoISet &  qoi_indices,
const libMesh::ParameterVector &  parameters_in,
libMesh::SensitivityData &  sensitivities 
) const
virtual

Reimplemented from GRINS::Solver.

Definition at line 164 of file steady_mesh_adaptive_solver.C.

References GRINS::SolverContext::system.

168  {
169  context.system->adjoint_qoi_parameter_sensitivity
170  (qoi_indices, parameters_in, sensitivities);
171  }
void GRINS::SteadyMeshAdaptiveSolver::check_qoi_error_option_consistency ( SolverContext context)
protected

Definition at line 65 of file steady_mesh_adaptive_solver.C.

References GRINS::MeshAdaptiveSolverBase::_error_estimator_options, GRINS::ErrorEstimatorOptions::compute_qoi_error_estimate(), and GRINS::SolverContext::error_estimator.

Referenced by solve().

66  {
67  // If we are computing QoI error estimates, we have to be using
68  // the Adjoint Refinement Error Estimator
70  if(context.error_estimator->type() != libMesh::ADJOINT_REFINEMENT)
71  {
72  std::string error_message = "You asked for QoI error estimates but did not use an Adjoint Refinement Error Estimator!\n";
73  error_message += "Please use the ADJOINT_REFINEMENT option for the estimator_type if you want QoI error estimates.\n";
74  std::cout<<error_message<<std::endl;
75  libmesh_error();
76  }
77  }
ErrorEstimatorOptions _error_estimator_options
void GRINS::SteadyMeshAdaptiveSolver::forward_qoi_parameter_sensitivity ( SolverContext context,
const libMesh::QoISet &  qoi_indices,
const libMesh::ParameterVector &  parameters_in,
libMesh::SensitivityData &  sensitivities 
) const
virtual

Reimplemented from GRINS::Solver.

Definition at line 174 of file steady_mesh_adaptive_solver.C.

References GRINS::SolverContext::equation_system, GRINS::SolverContext::output_residual_sensitivities, GRINS::SolverContext::output_solution_sensitivities, GRINS::SolverContext::system, and GRINS::SolverContext::vis.

178  {
179  context.system->forward_qoi_parameter_sensitivity
180  (qoi_indices, parameters_in, sensitivities);
181 
182  if( context.output_residual_sensitivities )
183  context.vis->output_residual_sensitivities
184  ( context.equation_system, context.system, parameters_in );
185 
186  if( context.output_solution_sensitivities )
187  context.vis->output_solution_sensitivities
188  ( context.equation_system, context.system, parameters_in );
189  }
void GRINS::SteadyMeshAdaptiveSolver::init_time_solver ( MultiphysicsSystem system)
protectedvirtual

Implements GRINS::Solver.

Definition at line 56 of file steady_mesh_adaptive_solver.C.

57  {
58  libMesh::SteadySolver* time_solver = new libMesh::SteadySolver( *(system) );
59 
60  system->time_solver = libMesh::UniquePtr<libMesh::TimeSolver>( time_solver );
61 
62  return;
63  }
void GRINS::SteadyMeshAdaptiveSolver::solve ( SolverContext context)
virtual
Todo:
This output cannot be toggled in the input file, but it should be able to be.

Implements GRINS::Solver.

Definition at line 79 of file steady_mesh_adaptive_solver.C.

References GRINS::MeshAdaptiveSolverBase::_error_estimator_options, GRINS::MeshAdaptiveSolverBase::_mesh_adaptivity_options, GRINS::MeshAdaptiveSolverBase::build_mesh_refinement(), GRINS::MeshAdaptiveSolverBase::check_for_convergence(), check_qoi_error_option_consistency(), GRINS::ErrorEstimatorOptions::compute_qoi_error_estimate(), GRINS::SolverContext::do_adjoint_solve, GRINS::SolverContext::equation_system, GRINS::SolverContext::error_estimator, GRINS::MeshAdaptiveSolverBase::estimate_error_for_amr(), GRINS::MeshAdaptivityOptions::max_refinement_steps(), GRINS::SolverContext::output_adjoint, GRINS::SolverContext::output_residual, GRINS::SolverContext::output_vis, GRINS::MeshAdaptiveSolverBase::perform_amr(), GRINS::SolverContext::postprocessing, GRINS::Solver::print_qoi(), GRINS::SolverContext::qoi_output, GRINS::Solver::steady_adjoint_solve(), GRINS::SolverContext::system, and GRINS::SolverContext::vis.

80  {
81  // Check that our error estimator options are consistent
83 
84  // Mesh and mesh refinement
85  libMesh::MeshBase& mesh = context.equation_system->get_mesh();
86  this->build_mesh_refinement( mesh );
87 
89  std::cout << "==========================================================" << std::endl
91  << " adaptive refinements" << std::endl
92  << "==========================================================" << std::endl;
93 
94  for ( unsigned int r_step = 0; r_step < _mesh_adaptivity_options.max_refinement_steps(); r_step++ )
95  {
96  std::cout << "==========================================================" << std::endl
97  << "Adaptive Refinement Step " << r_step << std::endl
98  << "==========================================================" << std::endl;
99 
100  // Solve the forward problem
101  context.system->solve();
102 
103  if( context.output_vis )
104  {
105  context.postprocessing->update_quantities( *(context.equation_system) );
106  context.vis->output( context.equation_system );
107  }
108 
109  // Solve adjoint system
110  if(context.do_adjoint_solve)
111  this->steady_adjoint_solve(context);
112 
113  if(context.output_adjoint)
114  context.vis->output_adjoint(context.equation_system, context.system);
115 
116  if( context.output_residual )
117  {
118  context.vis->output_residual( context.equation_system, context.system );
119  }
120 
121  // Now we construct the data structures for the mesh refinement process
122  libMesh::ErrorVector error;
123  this->estimate_error_for_amr( context, error );
124 
125  // Get the global error estimate if you can and are asked to
127  for(unsigned int i = 0; i != context.system->qoi.size(); i++)
128  {
129  libMesh::AdjointRefinementEstimator* adjoint_ref_error_estimator = libMesh::cast_ptr<libMesh::AdjointRefinementEstimator*>( context.error_estimator.get() );
130  std::cout<<"The error estimate for QoI("<<i<<") is: "<<adjoint_ref_error_estimator->get_global_QoI_error_estimate(i)<<std::endl;
131  }
132 
133  // Check for convergence of error
134  bool converged = this->check_for_convergence( context, error );
135 
136  if( converged )
137  {
138  // Break out of adaptive loop
139  std::cout << "==========================================================" << std::endl
140  << "Convergence detected!" << std::endl
141  << "==========================================================" << std::endl;
142  break;
143  }
144  else
145  {
146  // Only bother refining if we're on the last step.
148  {
149  this->perform_amr(context, error);
150 
151  // It's helpful to print the qoi along the way, but only do it if the user
152  // asks for it
153  if( context.qoi_output->output_qoi_set() )
154  this->print_qoi(context);
155  }
156  }
157 
158  } // r_step for-loop
159 
160  return;
161  }
void build_mesh_refinement(libMesh::MeshBase &mesh)
void check_qoi_error_option_consistency(SolverContext &context)
void perform_amr(SolverContext &context, const libMesh::ErrorVector &error)
void estimate_error_for_amr(SolverContext &context, libMesh::ErrorVector &error)
void steady_adjoint_solve(SolverContext &context)
Do steady version of adjoint solve.
Definition: solver.C:123
bool check_for_convergence(SolverContext &context, const libMesh::ErrorVector &error) const
ErrorEstimatorOptions _error_estimator_options
void print_qoi(SolverContext &context)
Definition: solver.C:155
unsigned int max_refinement_steps() const
MeshAdaptivityOptions _mesh_adaptivity_options

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