GRINS-0.8.0
mesh_adaptive_solver_base.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 // C++
26 #include <numeric>
27 #include <iomanip>
28 
29 // This class
31 #include "grins/solver_context.h"
32 
33 // libMesh
34 #include "libmesh/getpot.h"
35 #include "libmesh/mesh_base.h"
36 #include "libmesh/error_vector.h"
37 
38 namespace GRINS
39 {
41  : _error_estimator_options(input),
42  _mesh_adaptivity_options(input),
43  _refinement_type(INVALID)
44  {
46  }
47 
48  void MeshAdaptiveSolverBase::build_mesh_refinement( libMesh::MeshBase& mesh )
49  {
50  _mesh_refinement.reset( new libMesh::MeshRefinement( mesh ) );
60  _mesh_refinement->enforce_mismatch_limit_prior_to_refinement() = _mesh_adaptivity_options.enforce_mismatch_limit_prior_to_refinement();
62  }
63 
65  const MeshAdaptivityOptions& mesh_adaptivity_options,
67  {
68  const std::string& refinement_stategy = mesh_adaptivity_options.refinement_strategy();
69 
70  if( refinement_stategy == std::string("error_tolerance") )
71  {
72  refinement_type = ERROR_TOLERANCE;
73  }
74 
75  else if( refinement_stategy == std::string("nelem_target") )
76  {
77  refinement_type = N_ELEM_TARGET;
78 
79  if( !input.have_variable("MeshAdaptivity/nelem_target") )
80  {
81  std::cerr << "Error: Must specify nelem_target for" << std::endl
82  << " for error_tolerance refinement strategy."
83  << std::endl;
84 
85  libmesh_error();
86  }
87  }
88 
89  else if( refinement_stategy == std::string("error_fraction") )
90  {
91  refinement_type = ERROR_FRACTION;
92  }
93 
94  else if( refinement_stategy == std::string("elem_fraction") )
95  {
96  refinement_type = ELEM_FRACTION;
97  }
98 
99  else if( refinement_stategy == std::string("mean_std_dev") )
100  {
101  refinement_type = MEAN_STD_DEV;
102  }
103 
104  else
105  {
106  std::cerr << "Error: Invalid refinement strategy " << refinement_stategy << std::endl
107  << " Valid refinement strategy options are: absolute_global_tolerance" << std::endl
108  << " error_tolerance" << std::endl
109  << " nelem_target" << std::endl
110  << " error_fraction" << std::endl
111  << " elem_fraction" << std::endl
112  << " mean_std_dev" << std::endl;
113 
114  libmesh_error();
115  }
116 
117  return;
118  }
119 
121  const libMesh::ErrorVector& error ) const
122  {
123  std::cout << "==========================================================" << std::endl
124  << "Checking convergence" << std::endl
125  << "==========================================================" << std::endl;
126 
127  bool converged = false;
128 
129  libMesh::Real error_estimate = 0.0;
130 
131  if( context.do_adjoint_solve )
132  {
133  error_estimate = std::accumulate( error.begin(), error.end(), 0.0 );
134  }
135  else
136  {
137  error_estimate = error.l2_norm();
138  }
139 
140  std::cout << "==========================================================" << std::endl
141  << "Error estimate = " << error_estimate << std::endl
142  << "==========================================================" << std::endl;
143 
144  // For now, we just check the norm
145  if( std::fabs(error_estimate) <= _mesh_adaptivity_options.absolute_global_tolerance() )
146  {
147  converged = true;
148  }
149 
150  return converged;
151  }
152 
153  void MeshAdaptiveSolverBase::flag_elements_for_refinement( const libMesh::ErrorVector& error )
154  {
155  switch(_refinement_type)
156  {
157  case(ERROR_TOLERANCE):
158  {
159  _mesh_refinement->flag_elements_by_error_tolerance( error );
160  }
161  break;
162 
163  case(N_ELEM_TARGET):
164  {
165  _mesh_refinement->flag_elements_by_nelem_target( error );
166  }
167  break;
168 
169  case( ERROR_FRACTION ):
170  {
171  _mesh_refinement->flag_elements_by_error_fraction( error );
172  }
173  break;
174 
175  case( ELEM_FRACTION ):
176  {
177  _mesh_refinement->flag_elements_by_elem_fraction( error );
178  }
179  break;
180 
181  case( MEAN_STD_DEV ):
182  {
183  _mesh_refinement->flag_elements_by_mean_stddev( error );
184  }
185  break;
186 
187  case(INVALID):
188  {
189  std::cerr << "Error: Invalid refinement option!" << std::endl;
190  libmesh_error();
191  }
192  break;
193 
194  // Wat?!
195  default:
196  {
197  libmesh_error();
198  }
199 
200  } // switch(_refinement_type)
201 
202  return;
203  }
204 
205  void MeshAdaptiveSolverBase::estimate_error_for_amr( SolverContext& context, libMesh::ErrorVector& error )
206  {
207  std::cout << "==========================================================" << std::endl
208  << "Estimating error" << std::endl
209  << "==========================================================" << std::endl;
210  context.error_estimator->estimate_error( *context.system, error );
211 
212  libMesh::MeshBase& mesh = context.equation_system->get_mesh();
213 
214  // Plot error vector
216  {
217  error.plot_error( _mesh_adaptivity_options.error_plot_prefix()+".exo", mesh );
218  }
219  }
220 
221  void MeshAdaptiveSolverBase::perform_amr( SolverContext& context, const libMesh::ErrorVector& error )
222  {
223  libMesh::MeshBase& mesh = context.equation_system->get_mesh();
224 
225  std::cout << "==========================================================" << std::endl
226  << "Performing Mesh Refinement" << std::endl
227  << "==========================================================" << std::endl;
228 
229  this->flag_elements_for_refinement( error );
230  _mesh_refinement->refine_and_coarsen_elements();
231 
232  // Dont forget to reinit the system after each adaptive refinement!
233  context.equation_system->reinit();
234 
235  // This output cannot be toggled in the input file.
236  std::cout << "==========================================================" << std::endl
237  << "Refined mesh to " << std::setw(12) << mesh.n_active_elem()
238  << " active elements" << std::endl
239  << " " << std::setw(16) << context.system->n_active_dofs()
240  << " active dofs" << std::endl
241  << "==========================================================" << std::endl;
242  }
243 
244 } // end namespace GRINS
libMesh::Real coarsen_fraction() const
const std::string & refinement_strategy() const
void set_refinement_type(const GetPot &input, const MeshAdaptivityOptions &mesh_adaptivity_options, RefinementFlaggingType &refinement_type)
SharedPtr< libMesh::EquationSystems > equation_system
unsigned int edge_level_mismatch_limit() const
void build_mesh_refinement(libMesh::MeshBase &mesh)
libMesh::Real absolute_global_tolerance() const
unsigned int face_level_mismatch_limit() const
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
SharedPtr< libMesh::ErrorEstimator > error_estimator
GRINS namespace.
libMesh::Real coarsen_threshold() const
void flag_elements_for_refinement(const libMesh::ErrorVector &error)
libMesh::Real refine_fraction() const
Container for mesh adaptivity options.
GRINS::MultiphysicsSystem * system
bool enforce_mismatch_limit_prior_to_refinement() const
unsigned int node_level_mismatch_limit() const
Simple class to hold objects passed to Solver::solve.
libMesh::UniquePtr< libMesh::MeshRefinement > _mesh_refinement
const std::string & error_plot_prefix() const
MeshAdaptivityOptions _mesh_adaptivity_options

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