GRINS-0.8.0
steady_mesh_adaptive_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 // This class
27 
28 // C++
29 #include <iomanip>
30 
31 // GRINS
32 #include "grins/solver_context.h"
33 #include "grins/multiphysics_sys.h"
34 #include "grins/common.h"
35 
36 // libMesh
37 #include "libmesh/error_vector.h"
38 #include "libmesh/steady_solver.h"
39 #include "libmesh/adjoint_refinement_estimator.h"
40 
41 namespace GRINS
42 {
43 
45  : Solver(input),
46  MeshAdaptiveSolverBase( input )
47  {
48  return;
49  }
50 
52  {
53  return;
54  }
55 
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  }
64 
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  }
78 
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  }
162 
164  (SolverContext& context,
165  const libMesh::QoISet& qoi_indices,
166  const libMesh::ParameterVector& parameters_in,
167  libMesh::SensitivityData& sensitivities) const
168  {
169  context.system->adjoint_qoi_parameter_sensitivity
170  (qoi_indices, parameters_in, sensitivities);
171  }
172 
174  (SolverContext& context,
175  const libMesh::QoISet& qoi_indices,
176  const libMesh::ParameterVector& parameters_in,
177  libMesh::SensitivityData& sensitivities) const
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  }
190 
191 
192 } // end namespace GRINS
virtual void forward_qoi_parameter_sensitivity(SolverContext &context, const libMesh::QoISet &qoi_indices, const libMesh::ParameterVector &parameters_in, libMesh::SensitivityData &sensitivities) const
SharedPtr< libMesh::EquationSystems > equation_system
virtual void solve(SolverContext &context)
SharedPtr< PostProcessedQuantities< libMesh::Real > > postprocessing
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
SharedPtr< libMesh::ErrorEstimator > error_estimator
GRINS namespace.
SharedPtr< QoIOutput > qoi_output
GRINS::MultiphysicsSystem * system
SharedPtr< GRINS::Visualization > vis
Interface with libMesh for solving Multiphysics problems.
ErrorEstimatorOptions _error_estimator_options
Simple class to hold objects passed to Solver::solve.
void print_qoi(SolverContext &context)
Definition: solver.C:155
unsigned int max_refinement_steps() const
MeshAdaptivityOptions _mesh_adaptivity_options
virtual void adjoint_qoi_parameter_sensitivity(SolverContext &context, const libMesh::QoISet &qoi_indices, const libMesh::ParameterVector &parameters_in, libMesh::SensitivityData &sensitivities) const
virtual void init_time_solver(MultiphysicsSystem *system)

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