GRINS-0.8.0
steady_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 // This class
27 #include "grins/steady_solver.h"
28 
29 // GRINS
30 #include "grins/multiphysics_sys.h"
31 #include "grins/solver_context.h"
32 
33 // libMesh
34 #include "libmesh/auto_ptr.h"
35 #include "libmesh/dof_map.h"
36 #include "libmesh/getpot.h"
37 #include "libmesh/steady_solver.h"
38 #include "libmesh/linear_solver.h"
39 
40 namespace GRINS
41 {
42 
43  SteadySolver::SteadySolver( const GetPot& input )
44  : Solver( input )
45  {
46  return;
47  }
48 
50  {
51  return;
52  }
53 
55  {
56  libMesh::SteadySolver* time_solver = new libMesh::SteadySolver( *(system) );
57 
58  system->time_solver = libMesh::UniquePtr<libMesh::TimeSolver>(time_solver);
59  return;
60  }
61 
63  {
64  libmesh_assert( context.system );
65 
66  if( context.output_vis )
67  {
68  context.postprocessing->update_quantities( *(context.equation_system) );
69  context.vis->output( context.equation_system );
70  }
71 
72  context.system->solve();
73 
74  if ( context.print_scalars )
75  this->print_scalar_vars(context);
76 
77  if( context.do_adjoint_solve )
78  {
79  // Get the linear solver
80  libMesh::LinearSolver<libMesh::Number> *linear_solver = context.system->get_linear_solver();
81 
82  // Set ourselves to reuse the preconditioner
83  linear_solver->reuse_preconditioner(true);
84 
85  // Solve the adjoint problem
86  this->steady_adjoint_solve(context);
87 
88  // Go back to not reusing the preconditioner
89  linear_solver->reuse_preconditioner(false);
90  }
91 
92  if( context.output_adjoint )
93  context.vis->output_adjoint( context.equation_system, context.system );
94 
95  if( context.output_vis )
96  {
97  context.postprocessing->update_quantities( *(context.equation_system) );
98  context.vis->output( context.equation_system );
99  }
100 
101  if( context.output_residual ) context.vis->output_residual( context.equation_system, context.system );
102 
103  return;
104  }
105 
107  (SolverContext& context,
108  const libMesh::QoISet& qoi_indices,
109  const libMesh::ParameterVector& parameters_in,
110  libMesh::SensitivityData& sensitivities) const
111  {
112  // Get the linear solver
113  libMesh::LinearSolver<libMesh::Number> *linear_solver = context.system->get_linear_solver();
114 
115  // Set ourselves to reuse the preconditioner
116  linear_solver->reuse_preconditioner(true);
117 
118  context.system->adjoint_qoi_parameter_sensitivity
119  (qoi_indices, parameters_in, sensitivities);
120 
121  // Go back to not reusing the preconditioner
122  linear_solver->reuse_preconditioner(false);
123  }
124 
126  (SolverContext& context,
127  const libMesh::QoISet& qoi_indices,
128  const libMesh::ParameterVector& parameters_in,
129  libMesh::SensitivityData& sensitivities) const
130  {
131  context.system->forward_qoi_parameter_sensitivity
132  (qoi_indices, parameters_in, sensitivities);
133 
134  if( context.output_residual_sensitivities )
135  context.vis->output_residual_sensitivities
136  ( context.equation_system, context.system, parameters_in );
137 
138  if( context.output_solution_sensitivities )
139  context.vis->output_solution_sensitivities
140  ( context.equation_system, context.system, parameters_in );
141  }
142 
143 } // namespace GRINS
virtual void init_time_solver(GRINS::MultiphysicsSystem *system)
Definition: steady_solver.C:54
SharedPtr< libMesh::EquationSystems > equation_system
virtual void solve(SolverContext &context)
Definition: steady_solver.C:62
SharedPtr< PostProcessedQuantities< libMesh::Real > > postprocessing
void steady_adjoint_solve(SolverContext &context)
Do steady version of adjoint solve.
Definition: solver.C:123
GRINS namespace.
SteadySolver(const GetPot &input)
Definition: steady_solver.C:43
virtual void forward_qoi_parameter_sensitivity(SolverContext &context, const libMesh::QoISet &qoi_indices, const libMesh::ParameterVector &parameters_in, libMesh::SensitivityData &sensitivities) const
GRINS::MultiphysicsSystem * system
SharedPtr< GRINS::Visualization > vis
void print_scalar_vars(SolverContext &context)
Definition: solver.C:133
Interface with libMesh for solving Multiphysics problems.
virtual void adjoint_qoi_parameter_sensitivity(SolverContext &context, const libMesh::QoISet &qoi_indices, const libMesh::ParameterVector &parameters_in, libMesh::SensitivityData &sensitivities) const
Simple class to hold objects passed to Solver::solve.
virtual ~SteadySolver()
Definition: steady_solver.C:49

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