GRINS-0.6.0
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
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
 
virtual void initialize (const GetPot &input, std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system)
 
void steady_adjoint_solve (SolverContext &context)
 Do steady version of adjoint solve. More...
 

Protected Types

enum  RefinementFlaggingType {
  INVALID = 0, ERROR_TOLERANCE, N_ELEM_TARGET, ERROR_FRACTION,
  ELEM_FRACTION, MEAN_STD_DEV
}
 

Protected Member Functions

virtual void init_time_solver (MultiphysicsSystem *system)
 
void build_mesh_refinement (libMesh::MeshBase &mesh)
 
void set_refinement_type (const GetPot &input, RefinementFlaggingType &refinement_type)
 
bool check_for_convergence (SolverContext &context, const libMesh::ErrorVector &error) const
 
void flag_elements_for_refinement (const libMesh::ErrorVector &error)
 
void set_solver_options (libMesh::DiffSolver &solver)
 

Protected Attributes

unsigned int _max_refinement_steps
 
bool _coarsen_by_parents
 
libMesh::Real _absolute_global_tolerance
 
unsigned int _nelem_target
 
libMesh::Real _refine_fraction
 
libMesh::Real _coarsen_fraction
 
libMesh::Real _coarsen_threshold
 
bool _plot_cell_errors
 
std::string _error_plot_prefix
 
unsigned int _node_level_mismatch_limit
 
unsigned int _edge_level_mismatch_limit
 
unsigned int _face_level_mismatch_limit
 
bool _enforce_mismatch_limit_prior_to_refinement
 
RefinementFlaggingType _refinement_type
 
boost::scoped_ptr< libMesh::MeshRefinement > _mesh_refinement
 
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 steady_mesh_adaptive_solver.h.

Member Enumeration Documentation

Constructor & Destructor Documentation

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

Definition at line 40 of file steady_mesh_adaptive_solver.C.

41  : MeshAdaptiveSolverBase( input )
42  {
43  return;
44  }
GRINS::SteadyMeshAdaptiveSolver::~SteadyMeshAdaptiveSolver ( )
virtual

Definition at line 46 of file steady_mesh_adaptive_solver.C.

47  {
48  return;
49  }

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 168 of file steady_mesh_adaptive_solver.C.

References GRINS::SolverContext::system.

172  {
173  context.system->adjoint_qoi_parameter_sensitivity
174  (qoi_indices, parameters_in, sensitivities);
175  }
void GRINS::MeshAdaptiveSolverBase::build_mesh_refinement ( libMesh::MeshBase &  mesh)
protectedinherited

Definition at line 67 of file mesh_adaptive_solver_base.C.

References GRINS::MeshAdaptiveSolverBase::_absolute_global_tolerance, GRINS::MeshAdaptiveSolverBase::_coarsen_by_parents, GRINS::MeshAdaptiveSolverBase::_coarsen_fraction, GRINS::MeshAdaptiveSolverBase::_coarsen_threshold, GRINS::MeshAdaptiveSolverBase::_edge_level_mismatch_limit, GRINS::MeshAdaptiveSolverBase::_enforce_mismatch_limit_prior_to_refinement, GRINS::MeshAdaptiveSolverBase::_face_level_mismatch_limit, GRINS::MeshAdaptiveSolverBase::_mesh_refinement, GRINS::MeshAdaptiveSolverBase::_nelem_target, GRINS::MeshAdaptiveSolverBase::_node_level_mismatch_limit, and GRINS::MeshAdaptiveSolverBase::_refine_fraction.

Referenced by solve().

68  {
69  _mesh_refinement.reset( new libMesh::MeshRefinement( mesh ) );
70  _mesh_refinement->coarsen_by_parents() = _coarsen_by_parents;
71  _mesh_refinement->absolute_global_tolerance() = _absolute_global_tolerance;
72  _mesh_refinement->nelem_target() = _nelem_target;
73  _mesh_refinement->refine_fraction() = _refine_fraction;
74  _mesh_refinement->coarsen_fraction() = _coarsen_fraction;
75  _mesh_refinement->coarsen_threshold() = _coarsen_threshold;
76  _mesh_refinement->node_level_mismatch_limit() = _node_level_mismatch_limit;
77  _mesh_refinement->edge_level_mismatch_limit() = _edge_level_mismatch_limit;
78  _mesh_refinement->face_level_mismatch_limit() = _face_level_mismatch_limit;
79  _mesh_refinement->set_enforce_mismatch_limit_prior_to_refinement(_enforce_mismatch_limit_prior_to_refinement);
80 
81  return;
82  }
boost::scoped_ptr< libMesh::MeshRefinement > _mesh_refinement
bool GRINS::MeshAdaptiveSolverBase::check_for_convergence ( SolverContext context,
const libMesh::ErrorVector &  error 
) const
protectedinherited

Definition at line 148 of file mesh_adaptive_solver_base.C.

References GRINS::MeshAdaptiveSolverBase::_absolute_global_tolerance, and GRINS::SolverContext::do_adjoint_solve.

Referenced by solve().

150  {
151  bool converged = false;
152 
153  libMesh::Real error_estimate = 0.0;
154 
155  if( context.do_adjoint_solve )
156  {
157  error_estimate = std::accumulate( error.begin(), error.end(), 0.0 );
158  }
159  else
160  {
161  error_estimate = error.l2_norm();
162  }
163 
164  std::cout << "==========================================================" << std::endl
165  << "Error estimate = " << error_estimate << std::endl
166  << "==========================================================" << std::endl;
167 
168  // For now, we just check the norm
169  if( std::fabs(error_estimate) <= _absolute_global_tolerance )
170  {
171  converged = true;
172  }
173 
174  return converged;
175  }
void GRINS::MeshAdaptiveSolverBase::flag_elements_for_refinement ( const libMesh::ErrorVector &  error)
protectedinherited

Definition at line 177 of file mesh_adaptive_solver_base.C.

References GRINS::MeshAdaptiveSolverBase::_mesh_refinement, GRINS::MeshAdaptiveSolverBase::_refinement_type, GRINS::MeshAdaptiveSolverBase::ELEM_FRACTION, GRINS::MeshAdaptiveSolverBase::ERROR_FRACTION, GRINS::MeshAdaptiveSolverBase::ERROR_TOLERANCE, GRINS::MeshAdaptiveSolverBase::INVALID, GRINS::MeshAdaptiveSolverBase::MEAN_STD_DEV, and GRINS::MeshAdaptiveSolverBase::N_ELEM_TARGET.

Referenced by solve().

178  {
179  switch(_refinement_type)
180  {
181  case(ERROR_TOLERANCE):
182  {
183  _mesh_refinement->flag_elements_by_error_tolerance( error );
184  }
185  break;
186 
187  case(N_ELEM_TARGET):
188  {
189  _mesh_refinement->flag_elements_by_nelem_target( error );
190  }
191  break;
192 
193  case( ERROR_FRACTION ):
194  {
195  _mesh_refinement->flag_elements_by_error_fraction( error );
196  }
197  break;
198 
199  case( ELEM_FRACTION ):
200  {
201  _mesh_refinement->flag_elements_by_elem_fraction( error );
202  }
203  break;
204 
205  case( MEAN_STD_DEV ):
206  {
207  _mesh_refinement->flag_elements_by_mean_stddev( error );
208  }
209  break;
210 
211  case(INVALID):
212  {
213  std::cerr << "Error: Invalid refinement option!" << std::endl;
214  libmesh_error();
215  }
216  break;
217 
218  // Wat?!
219  default:
220  {
221  libmesh_error();
222  }
223 
224  } // switch(_refinement_type)
225 
226  return;
227  }
boost::scoped_ptr< libMesh::MeshRefinement > _mesh_refinement
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 178 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.

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

Implements GRINS::Solver.

Definition at line 51 of file steady_mesh_adaptive_solver.C.

52  {
53  libMesh::SteadySolver* time_solver = new libMesh::SteadySolver( *(system) );
54 
55  system->time_solver = libMesh::AutoPtr<libMesh::TimeSolver>( time_solver );
56 
57  return;
58  }
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::MeshAdaptiveSolverBase::set_refinement_type ( const GetPot &  input,
MeshAdaptiveSolverBase::RefinementFlaggingType refinement_type 
)
protectedinherited

Definition at line 84 of file mesh_adaptive_solver_base.C.

References GRINS::MeshAdaptiveSolverBase::ELEM_FRACTION, GRINS::MeshAdaptiveSolverBase::ERROR_FRACTION, GRINS::MeshAdaptiveSolverBase::ERROR_TOLERANCE, GRINS::MeshAdaptiveSolverBase::MEAN_STD_DEV, and GRINS::MeshAdaptiveSolverBase::N_ELEM_TARGET.

Referenced by GRINS::MeshAdaptiveSolverBase::MeshAdaptiveSolverBase().

86  {
87  std::string refinement_stategy = input("MeshAdaptivity/refinement_strategy", "elem_fraction" );
88 
89  if( !input.have_variable("MeshAdaptivity/absolute_global_tolerance" ) )
90  {
91  std::cerr << "Error: Must specify absolute_global_tolerance for" << std::endl
92  << " adaptive refinement algorithms."
93  << std::endl;
94 
95  libmesh_error();
96  }
97 
98  if( refinement_stategy == std::string("error_tolerance") )
99  {
100  refinement_type = ERROR_TOLERANCE;
101  }
102 
103  else if( refinement_stategy == std::string("nelem_target") )
104  {
105  refinement_type = N_ELEM_TARGET;
106 
107  if( !input.have_variable("MeshAdaptivity/nelem_target") )
108  {
109  std::cerr << "Error: Must specify nelem_target for" << std::endl
110  << " for error_tolerance refinement strategy."
111  << std::endl;
112 
113  libmesh_error();
114  }
115  }
116 
117  else if( refinement_stategy == std::string("error_fraction") )
118  {
119  refinement_type = ERROR_FRACTION;
120  }
121 
122  else if( refinement_stategy == std::string("elem_fraction") )
123  {
124  refinement_type = ELEM_FRACTION;
125  }
126 
127  else if( refinement_stategy == std::string("mean_std_dev") )
128  {
129  refinement_type = MEAN_STD_DEV;
130  }
131 
132  else
133  {
134  std::cerr << "Error: Invalid refinement strategy " << refinement_stategy << std::endl
135  << " Valid refinement strategy options are: absolute_global_tolerance" << std::endl
136  << " error_tolerance" << std::endl
137  << " nelem_target" << std::endl
138  << " error_fraction" << std::endl
139  << " elem_fraction" << std::endl
140  << " mean_std_dev" << std::endl;
141 
142  libmesh_error();
143  }
144 
145  return;
146  }
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::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 60 of file steady_mesh_adaptive_solver.C.

References GRINS::MeshAdaptiveSolverBase::_error_plot_prefix, GRINS::MeshAdaptiveSolverBase::_max_refinement_steps, GRINS::MeshAdaptiveSolverBase::_mesh_refinement, GRINS::MeshAdaptiveSolverBase::_plot_cell_errors, GRINS::MeshAdaptiveSolverBase::build_mesh_refinement(), GRINS::MeshAdaptiveSolverBase::check_for_convergence(), GRINS::SolverContext::do_adjoint_solve, GRINS::SolverContext::equation_system, GRINS::SolverContext::error_estimator, GRINS::MeshAdaptiveSolverBase::flag_elements_for_refinement(), GRINS::SolverContext::output_adjoint, GRINS::CompositeQoI::output_qoi(), GRINS::SolverContext::output_residual, GRINS::SolverContext::output_vis, GRINS::SolverContext::postprocessing, GRINS::SolverContext::print_qoi, GRINS::Solver::steady_adjoint_solve(), GRINS::SolverContext::system, and GRINS::SolverContext::vis.

61  {
62  // Mesh and mesh refinement
63  libMesh::MeshBase& mesh = context.equation_system->get_mesh();
64  this->build_mesh_refinement( mesh );
65 
67  std::cout << "==========================================================" << std::endl
68  << "Performing " << this->_max_refinement_steps << " adaptive refinements" << std::endl
69  << "==========================================================" << std::endl;
70 
71  // GRVY timers contained in here (if enabled)
72  for ( unsigned int r_step = 0; r_step < this->_max_refinement_steps; r_step++ )
73  {
74  std::cout << "==========================================================" << std::endl
75  << "Adaptive Refinement Step " << r_step << std::endl
76  << "==========================================================" << std::endl;
77 
78  // Solve the forward problem
79  context.system->solve();
80 
81  if( context.output_vis )
82  {
83  context.postprocessing->update_quantities( *(context.equation_system) );
84  context.vis->output( context.equation_system );
85  }
86 
87  // Solve adjoint system
88  if(context.do_adjoint_solve)
89  this->steady_adjoint_solve(context);
90 
91  if(context.output_adjoint)
92  context.vis->output_adjoint(context.equation_system, context.system);
93 
94  if( context.output_residual )
95  {
96  context.vis->output_residual( context.equation_system, context.system );
97  }
98 
99  // Now we construct the data structures for the mesh refinement process
100  libMesh::ErrorVector error;
101 
102  std::cout << "==========================================================" << std::endl
103  << "Estimating error" << std::endl
104  << "==========================================================" << std::endl;
105  context.error_estimator->estimate_error( *context.system, error );
106 
107  // Plot error vector
108  if( this->_plot_cell_errors )
109  {
110  error.plot_error( this->_error_plot_prefix+".exo", mesh );
111  }
112 
113  // Check for convergence of error
114  std::cout << "==========================================================" << std::endl
115  << "Checking convergence" << std::endl
116  << "==========================================================" << std::endl;
117  bool converged = this->check_for_convergence( context, error );
118 
119  if( converged )
120  {
121  // Break out of adaptive loop
122  std::cout << "==========================================================" << std::endl
123  << "Convergence detected!" << std::endl
124  << "==========================================================" << std::endl;
125  break;
126  }
127  else
128  {
129  // Only bother refining if we're on the last step.
130  if( r_step < this->_max_refinement_steps -1 )
131  {
132  std::cout << "==========================================================" << std::endl
133  << "Performing Mesh Refinement" << std::endl
134  << "==========================================================" << std::endl;
135 
136  this->flag_elements_for_refinement( error );
137  _mesh_refinement->refine_and_coarsen_elements();
138 
139  // Dont forget to reinit the system after each adaptive refinement!
140  context.equation_system->reinit();
141 
142  // This output cannot be toggled in the input file.
143  std::cout << "==========================================================" << std::endl
144  << "Refined mesh to " << std::setw(12) << mesh.n_active_elem()
145  << " active elements" << std::endl
146  << " " << std::setw(16) << context.system->n_active_dofs()
147  << " active dofs" << std::endl
148  << "==========================================================" << std::endl;
149 
150  // It's helpful to print the qoi along the way, but only do it if the user
151  // asks for it
152  if( context.print_qoi )
153  {
154  context.system->assemble_qoi();
155  const CompositeQoI* my_qoi = libMesh::libmesh_cast_ptr<const CompositeQoI*>(context.system->get_qoi());
156  my_qoi->output_qoi( std::cout );
157  std::cout << std::endl;
158  }
159  }
160  }
161 
162  } // r_step for-loop
163 
164  return;
165  }
boost::scoped_ptr< libMesh::MeshRefinement > _mesh_refinement
void build_mesh_refinement(libMesh::MeshBase &mesh)
void steady_adjoint_solve(SolverContext &context)
Do steady version of adjoint solve.
Definition: grins_solver.C:105
bool check_for_convergence(SolverContext &context, const libMesh::ErrorVector &error) const
void flag_elements_for_refinement(const libMesh::ErrorVector &error)
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 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

libMesh::Real GRINS::MeshAdaptiveSolverBase::_absolute_global_tolerance
protectedinherited
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().

bool GRINS::MeshAdaptiveSolverBase::_coarsen_by_parents
protectedinherited
libMesh::Real GRINS::MeshAdaptiveSolverBase::_coarsen_fraction
protectedinherited
libMesh::Real GRINS::MeshAdaptiveSolverBase::_coarsen_threshold
protectedinherited
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().

unsigned int GRINS::MeshAdaptiveSolverBase::_edge_level_mismatch_limit
protectedinherited
bool GRINS::MeshAdaptiveSolverBase::_enforce_mismatch_limit_prior_to_refinement
protectedinherited
std::string GRINS::MeshAdaptiveSolverBase::_error_plot_prefix
protectedinherited

Definition at line 76 of file mesh_adaptive_solver_base.h.

Referenced by solve().

unsigned int GRINS::MeshAdaptiveSolverBase::_face_level_mismatch_limit
protectedinherited
double GRINS::Solver::_initial_linear_tolerance
protectedinherited

Definition at line 104 of file grins_solver.h.

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

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().

unsigned int GRINS::MeshAdaptiveSolverBase::_max_refinement_steps
protectedinherited

Definition at line 68 of file mesh_adaptive_solver_base.h.

Referenced by solve().

boost::scoped_ptr<libMesh::MeshRefinement> GRINS::MeshAdaptiveSolverBase::_mesh_refinement
protectedinherited
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::MeshAdaptiveSolverBase::_nelem_target
protectedinherited
std::map< std::string, GRINS::NBCContainer > GRINS::Solver::_neumann_bc_funcs
protectedinherited

Definition at line 117 of file grins_solver.h.

unsigned int GRINS::MeshAdaptiveSolverBase::_node_level_mismatch_limit
protectedinherited
bool GRINS::MeshAdaptiveSolverBase::_plot_cell_errors
protectedinherited

Definition at line 75 of file mesh_adaptive_solver_base.h.

Referenced by solve().

libMesh::Real GRINS::MeshAdaptiveSolverBase::_refine_fraction
protectedinherited
RefinementFlaggingType GRINS::MeshAdaptiveSolverBase::_refinement_type
protectedinherited
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().


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