GRINS-0.6.0
grins_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-2015 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
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 
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::AutoPtr<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  // GRVY timers contained in here (if enabled)
73  context.system->solve();
74 
75  if ( context.print_scalars )
76  for (unsigned int v=0; v != context.system->n_vars(); ++v)
77  if (context.system->variable(v).type().family ==
78  libMesh::SCALAR)
79  {
80  std::cout << context.system->variable_name(v) <<
81  " = {";
82  std::vector<libMesh::dof_id_type> scalar_indices;
83  context.system->get_dof_map().SCALAR_dof_indices
84  (scalar_indices, v);
85  if (scalar_indices.size())
86  std::cout <<
87  context.system->current_solution(scalar_indices[0]);
88  for (unsigned int i=1; i < scalar_indices.size();
89  ++i)
90  std::cout << ", " <<
91  context.system->current_solution(scalar_indices[i]);
92  std::cout << '}' << std::endl;
93  }
94 
95  if( context.do_adjoint_solve )
96  this->steady_adjoint_solve(context);
97 
98  if( context.output_adjoint )
99  context.vis->output_adjoint( context.equation_system, context.system );
100 
101  if( context.output_vis )
102  {
103  context.postprocessing->update_quantities( *(context.equation_system) );
104  context.vis->output( context.equation_system );
105  }
106 
107  if( context.output_residual ) context.vis->output_residual( context.equation_system, context.system );
108 
109  return;
110  }
111 
113  (SolverContext& context,
114  const libMesh::QoISet& qoi_indices,
115  const libMesh::ParameterVector& parameters_in,
116  libMesh::SensitivityData& sensitivities) const
117  {
118  context.system->adjoint_qoi_parameter_sensitivity
119  (qoi_indices, parameters_in, sensitivities);
120  }
121 
123  (SolverContext& context,
124  const libMesh::QoISet& qoi_indices,
125  const libMesh::ParameterVector& parameters_in,
126  libMesh::SensitivityData& sensitivities) const
127  {
128  context.system->forward_qoi_parameter_sensitivity
129  (qoi_indices, parameters_in, sensitivities);
130 
131  if( context.output_residual_sensitivities )
132  context.vis->output_residual_sensitivities
133  ( context.equation_system, context.system, parameters_in );
134 
135  if( context.output_solution_sensitivities )
136  context.vis->output_solution_sensitivities
137  ( context.equation_system, context.system, parameters_in );
138  }
139 
140 } // namespace GRINS
virtual void init_time_solver(GRINS::MultiphysicsSystem *system)
virtual void solve(SolverContext &context)
void steady_adjoint_solve(SolverContext &context)
Do steady version of adjoint solve.
Definition: grins_solver.C:105
std::tr1::shared_ptr< libMesh::EquationSystems > equation_system
GRINS namespace.
SteadySolver(const GetPot &input)
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
Interface with libMesh for solving Multiphysics problems.
std::tr1::shared_ptr< PostProcessedQuantities< libMesh::Real > > postprocessing
std::tr1::shared_ptr< GRINS::Visualization > vis
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.

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1