GRINS-0.8.0
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 
26 // C++
27 #include <iostream>
28 #include <iomanip>
29 
30 // This class
31 #include "grins/solver.h"
32 
33 // GRINS
34 #include "grins/multiphysics_sys.h"
35 #include "grins/solver_context.h"
36 #include "grins/composite_qoi.h"
37 
38 // libMesh
39 #include "libmesh/getpot.h"
40 #include "libmesh/fem_system.h"
41 #include "libmesh/diff_solver.h"
42 #include "libmesh/newton_solver.h"
43 #include "libmesh/dof_map.h"
44 
45 namespace GRINS
46 {
47 
48  Solver::Solver( const GetPot& input )
49  : _max_nonlinear_iterations( input("linear-nonlinear-solver/max_nonlinear_iterations", 10 ) ),
50  _relative_step_tolerance( input("linear-nonlinear-solver/relative_step_tolerance", 1.e-6 ) ),
51  _absolute_step_tolerance( input("linear-nonlinear-solver/absolute_step_tolerance", 0.0 ) ),
52  _relative_residual_tolerance( input("linear-nonlinear-solver/relative_residual_tolerance", 1.e-15 ) ),
53  _absolute_residual_tolerance( input("linear-nonlinear-solver/absolute_residual_tolerance", 0.0 ) ),
54  _initial_linear_tolerance( input("linear-nonlinear-solver/initial_linear_tolerance", 1.e-3 ) ),
55  _minimum_linear_tolerance( input("linear-nonlinear-solver/minimum_linear_tolerance", 1.e-3 ) ),
56  _max_linear_iterations( input("linear-nonlinear-solver/max_linear_iterations", 500 ) ),
57  _continue_after_backtrack_failure( input("linear-nonlinear-solver/continue_after_backtrack_failure", false ) ),
58  _continue_after_max_iterations( input("linear-nonlinear-solver/continue_after_max_iterations", false ) ),
59  _require_residual_reduction( input("linear-nonlinear-solver/require_residual_reduction", true ) ),
60  _solver_quiet( input("screen-options/solver_quiet", false ) ),
61  _solver_verbose( input("screen-options/solver_verbose", false ) )
62  {
63  return;
64  }
65 
66 
68  {
69  return;
70  }
71 
72  void Solver::initialize( const GetPot& /*input*/,
73  SharedPtr<libMesh::EquationSystems> equation_system,
74  MultiphysicsSystem* system )
75  {
76 
77  // Defined in subclasses depending on the solver used.
78  this->init_time_solver(system);
79 
80  // Initialize the system
81  equation_system->init();
82 
83  // Get diff solver to set options
84  libMesh::DiffSolver &solver = *(system->time_solver->diff_solver().get());
85 
86  // Set linear/nonlinear solver options
87  this->set_solver_options( solver );
88 
89  return;
90  }
91 
92  void Solver::set_solver_options( libMesh::DiffSolver& solver )
93  {
94  solver.quiet = this->_solver_quiet;
95  solver.verbose = this->_solver_verbose;
96  solver.max_nonlinear_iterations = this->_max_nonlinear_iterations;
97  solver.relative_step_tolerance = this->_relative_step_tolerance;
98  solver.absolute_step_tolerance = this->_absolute_step_tolerance;
99  solver.relative_residual_tolerance = this->_relative_residual_tolerance;
100  solver.absolute_residual_tolerance = this->_absolute_residual_tolerance;
101  solver.max_linear_iterations = this->_max_linear_iterations;
102  solver.continue_after_backtrack_failure = this->_continue_after_backtrack_failure;
103  solver.initial_linear_tolerance = this->_initial_linear_tolerance;
104  solver.minimum_linear_tolerance = this->_minimum_linear_tolerance;
105  solver.continue_after_max_iterations = this->_continue_after_max_iterations;
106  if(dynamic_cast<libMesh::NewtonSolver*>(&solver))
107  {
108  dynamic_cast<libMesh::NewtonSolver&>(solver).require_residual_reduction = this->_require_residual_reduction;
109  }
110  else
111  {
112  // If the user tried to set require_residual_reduction flag to false
113  // despite not having a NewtonSolver spit out a warning
114  if(this->_require_residual_reduction == false)
115  {
116  libmesh_warning("GRINS can't change require_residual_reduction when not using NewtonSolver!");
117  }
118  }
119 
120  return;
121  }
122 
124  {
125  libMesh::out << "==========================================================" << std::endl
126  << "Solving adjoint problem." << std::endl
127  << "==========================================================" << std::endl;
128 
129  context.system->adjoint_solve();
130  context.system->set_adjoint_already_solved(true);
131  }
132 
134  {
135  for (unsigned int v=0; v != context.system->n_vars(); ++v)
136  if (context.system->variable(v).type().family ==
137  libMesh::SCALAR)
138  {
139  std::cout << context.system->variable_name(v) <<
140  " = {";
141  std::vector<libMesh::dof_id_type> scalar_indices;
142  context.system->get_dof_map().SCALAR_dof_indices
143  (scalar_indices, v);
144  if (scalar_indices.size())
145  std::cout <<
146  context.system->current_solution(scalar_indices[0]);
147  for (unsigned int i=1; i < scalar_indices.size();
148  ++i)
149  std::cout << ", " <<
150  context.system->current_solution(scalar_indices[i]);
151  std::cout << '}' << std::endl;
152  }
153  }
154 
156  {
157  context.system->assemble_qoi();
158  const CompositeQoI* my_qoi = libMesh::cast_ptr<const CompositeQoI*>(context.system->get_qoi());
159  context.qoi_output->output_qois(*my_qoi, context.system->comm());
160  }
161 
162 } // namespace GRINS
double _absolute_residual_tolerance
Definition: solver.h:100
bool _solver_verbose
Definition: solver.h:110
virtual void init_time_solver(GRINS::MultiphysicsSystem *system)=0
double _relative_residual_tolerance
Definition: solver.h:98
Solver(const GetPot &input)
Definition: solver.C:48
unsigned int _max_linear_iterations
Definition: solver.h:103
void set_solver_options(libMesh::DiffSolver &solver)
Definition: solver.C:92
double _relative_step_tolerance
Definition: solver.h:92
void steady_adjoint_solve(SolverContext &context)
Do steady version of adjoint solve.
Definition: solver.C:123
bool _continue_after_backtrack_failure
Definition: solver.h:104
GRINS namespace.
double _minimum_linear_tolerance
Definition: solver.h:102
SharedPtr< QoIOutput > qoi_output
virtual void initialize(const GetPot &input, SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system)
Definition: solver.C:72
double _absolute_step_tolerance
Definition: solver.h:93
GRINS::MultiphysicsSystem * system
void print_scalar_vars(SolverContext &context)
Definition: solver.C:133
Interface with libMesh for solving Multiphysics problems.
bool _continue_after_max_iterations
Definition: solver.h:105
virtual ~Solver()
Definition: solver.C:67
Simple class to hold objects passed to Solver::solve.
void print_qoi(SolverContext &context)
Definition: solver.C:155
unsigned int _max_nonlinear_iterations
Definition: solver.h:91
bool _require_residual_reduction
Definition: solver.h:106
double _initial_linear_tolerance
Definition: solver.h:101
bool _solver_quiet
Definition: solver.h:109

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