GRINS-0.7.0
unsteady_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-2016 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 #include "libmesh/steady_solver.h"
37 
38 namespace GRINS
39 {
40 
42  ( const GetPot& input,
43  const libMesh::Parallel::Communicator &comm )
44  : Visualization(input, comm)
45  {
46  return;
47  }
48 
50  {
51  return;
52  }
53 
55  ( SharedPtr<libMesh::EquationSystems> equation_system,
56  MultiphysicsSystem* system,
57  const unsigned int time_step,
58  const libMesh::Real time )
59  {
60  std::stringstream suffix;
61  suffix << time_step;
62 
63  std::string filename = this->_vis_output_file_prefix+"_unsteady_residual";
64 
65  filename+="."+suffix.str();
66 
67  // For the unsteady residual, we just want to evaluate F(u) from
68  // dU/dt = F(u). What we do is swap out the time solver to a
69  // SteadySolver and reassemble the residual. Then, we'll need to swap
70  // the solution and the rhs vector stashed in the system. Once we're done,
71  // we'll reset the time solver pointer back to the original guy.
72 
73  libMesh::TimeSolver* prev_time_solver = system->time_solver.get();
74 
75  libMesh::SteadySolver* steady_solver = new libMesh::SteadySolver( *(system) );
76 
77  system->time_solver = libMesh::UniquePtr<libMesh::TimeSolver>(steady_solver);
78 
79  system->assembly( true /*residual*/, false /*jacobian*/ );
80  system->rhs->close();
81 
82  // Swap solution with newly computed residual
83  system->solution->swap( *(system->rhs) );
84  // Update equation systems
85  equation_system->update();
86 
87  this->dump_visualization( equation_system, filename, time );
88 
89  // Now swap back and reupdate
90  system->solution->swap( *(system->rhs) );
91  equation_system->update();
92 
93  system->time_solver.reset(prev_time_solver);
94  }
95 
97  (SharedPtr<libMesh::EquationSystems> equation_system,
98  MultiphysicsSystem* system,
99  const libMesh::ParameterVector & params,
100  const unsigned int time_step,
101  const libMesh::Real time )
102  {
103  for (unsigned int p=0; p != params.size(); ++p)
104  {
105  std::stringstream suffix;
106  suffix << time_step;
107 
108  std::stringstream pstr;
109  pstr << p;
110 
111  std::string filename =
112  this->_vis_output_file_prefix + "_unsteady_dRdp" +
113  pstr.str() + '.' + suffix.str();
114 
115  // Swap solution with precomputed sensitivity rhs
116  system->solution->swap(system->get_sensitivity_rhs(p));
117  equation_system->update();
118 
119  this->dump_visualization( equation_system, filename, time );
120 
121  // Now swap back and reupdate
122  system->solution->swap(system->get_sensitivity_rhs(p));
123  equation_system->update();
124  }
125  }
126 
128  ( SharedPtr<libMesh::EquationSystems> equation_system,
129  MultiphysicsSystem* system,
130  const unsigned int time_step,
131  const libMesh::Real time )
132  {
133  std::stringstream suffix;
134  suffix << time_step;
135 
136  const libMesh::DifferentiableQoI* raw_qoi = system->get_qoi();
137  const CompositeQoI* qoi = dynamic_cast<const CompositeQoI*>( raw_qoi );
138 
139  unsigned int n_qois = qoi->n_qois();
140 
141  for( unsigned int q = 0; q < n_qois; q++ )
142  {
143  libMesh::NumericVector<libMesh::Number>& dual_solution = system->get_adjoint_solution(q);
144 
145  const std::string& qoi_name = qoi->get_qoi(q).name();
146  std::string filename = this->_vis_output_file_prefix+"_unsteady_adjoint_"+qoi_name;
147  filename+="."+suffix.str();
148 
149  system->solution->swap( dual_solution );
150  equation_system->update();
151 
152  this->dump_visualization( equation_system, filename, time );
153 
154  // Now swap back and reupdate
155  system->solution->swap( dual_solution );
156  equation_system->update();
157  }
158  }
159 
161  (SharedPtr<libMesh::EquationSystems> equation_system,
162  MultiphysicsSystem* system,
163  const libMesh::ParameterVector & params,
164  const unsigned int time_step,
165  const libMesh::Real time )
166  {
167  for (unsigned int p=0; p != params.size(); ++p)
168  {
169  std::stringstream suffix;
170  suffix << time_step;
171 
172  std::stringstream pstr;
173  pstr << p;
174 
175  std::string filename =
176  this->_vis_output_file_prefix + "_unsteady_dudp" +
177  pstr.str() + '.' + suffix.str();
178 
179  // Swap solution with precomputed sensitivity solution
180  system->solution->swap(system->get_sensitivity_solution(p));
181  equation_system->update();
182 
183  this->dump_visualization( equation_system, filename, time );
184 
185  // Now swap back and reupdate
186  system->solution->swap(system->get_sensitivity_solution(p));
187  equation_system->update();
188  }
189  }
190 
191 
192 } // namespace GRINS
const std::string & name() const
Returns the name of this QoI.
Definition: qoi_base.h:139
virtual void output_residual(SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const unsigned int time_step, const libMesh::Real time)
virtual void output_residual_sensitivities(SharedPtr< 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
virtual void output_solution_sensitivities(SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const libMesh::ParameterVector &params, const unsigned int time_step, const libMesh::Real time)
GRINS namespace.
UnsteadyVisualization(const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
virtual void output_adjoint(SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const unsigned int time_step, const libMesh::Real time)
Interface with libMesh for solving Multiphysics problems.
unsigned int n_qois() const

Generated on Thu Jun 2 2016 21:52:28 for GRINS-0.7.0 by  doxygen 1.8.10