GRINS-0.6.0
steady_visualization.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/composite_qoi.h"
32 
33 // libMesh
34 #include "libmesh/getpot.h"
35 #include "libmesh/parameter_vector.h"
36 
37 namespace GRINS
38 {
39 
41  ( const GetPot& input,
42  const libMesh::Parallel::Communicator &comm )
43  : Visualization(input, comm)
44  {
45  return;
46  }
47 
49  {
50  return;
51  }
52 
53  void SteadyVisualization::output_residual( std::tr1::shared_ptr<libMesh::EquationSystems> equation_system,
54  MultiphysicsSystem* system,
55  const unsigned int,
56  const libMesh::Real )
57  {
58  std::string filename = this->_vis_output_file_prefix+"_residual";
59 
60  // Idea is that this->rhs stashes the residual. Thus, when we swap
61  // with the solution, we should be dumping the residual. Then, we swap
62  // back once we're done outputting.
63 
64  // Swap solution with computed residual
65  system->solution->swap( *(system->rhs) );
66  equation_system->update();
67 
68  this->dump_visualization( equation_system, filename, 0.0 );
69 
70  // Now swap back and reupdate
71  system->solution->swap( *(system->rhs) );
72  equation_system->update();
73 
74  return;
75  }
76 
78  (std::tr1::shared_ptr<libMesh::EquationSystems> equation_system,
79  MultiphysicsSystem* system,
80  const libMesh::ParameterVector & params,
81  const unsigned int,
82  const libMesh::Real )
83  {
84  for (unsigned int p=0; p != params.size(); ++p)
85  {
86  std::stringstream pstr;
87  pstr << p;
88 
89  std::string filename =
90  this->_vis_output_file_prefix+"_dRdp"+pstr.str();
91 
92  // Swap solution with precomputed sensitivity rhs
93  system->solution->swap(system->get_sensitivity_rhs(p));
94  equation_system->update();
95 
96  this->dump_visualization( equation_system, filename, 0.0 );
97 
98  // Now swap back and reupdate
99  system->solution->swap(system->get_sensitivity_rhs(p));
100  equation_system->update();
101  }
102  }
103 
104  void SteadyVisualization::output_adjoint( std::tr1::shared_ptr<libMesh::EquationSystems> equation_system,
105  MultiphysicsSystem* system,
106  const unsigned int /*time_step*/,
107  const libMesh::Real /*time*/ )
108  {
109  const libMesh::DifferentiableQoI* raw_qoi = system->get_qoi();
110  const CompositeQoI* qoi = dynamic_cast<const CompositeQoI*>( raw_qoi );
111 
112  unsigned int n_qois = qoi->n_qois();
113 
114  for( unsigned int q = 0; q < n_qois; q++ )
115  {
116  libMesh::NumericVector<libMesh::Number>& dual_solution = system->get_adjoint_solution(q);
117 
118  const std::string& qoi_name = qoi->get_qoi(q).name();
119  std::string filename = this->_vis_output_file_prefix+"_adjoint_"+qoi_name;
120 
121  system->solution->swap( dual_solution );
122  equation_system->update();
123 
124  this->dump_visualization( equation_system, filename, 0.0 );
125 
126  // Now swap back and reupdate
127  system->solution->swap( dual_solution );
128  equation_system->update();
129  }
130  }
131 
133  (std::tr1::shared_ptr<libMesh::EquationSystems> equation_system,
134  MultiphysicsSystem* system,
135  const libMesh::ParameterVector & params,
136  const unsigned int,
137  const libMesh::Real )
138  {
139  for (unsigned int p=0; p != params.size(); ++p)
140  {
141  std::stringstream pstr;
142  pstr << p;
143 
144  std::string filename =
145  this->_vis_output_file_prefix+"_dudp"+pstr.str();
146 
147  // Swap solution with precomputed sensitivity solution
148  system->solution->swap(system->get_sensitivity_solution(p));
149  equation_system->update();
150 
151  this->dump_visualization( equation_system, filename, 0.0 );
152 
153  // Now swap back and reupdate
154  system->solution->swap(system->get_sensitivity_solution(p));
155  equation_system->update();
156  }
157  }
158 
159 
160 } // namespace GRINS
virtual void output_residual(std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const unsigned int time_step, const libMesh::Real time)
void dump_visualization(std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, const std::string &filename_prefix, const libMesh::Real time)
const std::string & name() const
Returns the name of this QoI.
Definition: qoi_base.h:143
virtual void output_residual_sensitivities(std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const libMesh::ParameterVector &params, const unsigned int time_step, const libMesh::Real time)
const QoIBase & get_qoi(unsigned int qoi_index) const
std::string _vis_output_file_prefix
GRINS namespace.
virtual void output_solution_sensitivities(std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const libMesh::ParameterVector &params, const unsigned int time_step, const libMesh::Real time)
Interface with libMesh for solving Multiphysics problems.
SteadyVisualization(const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
virtual void output_adjoint(std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const unsigned int time_step, const libMesh::Real time)
unsigned int n_qois() const

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