GRINS-0.6.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-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 #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  ( std::tr1::shared_ptr<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::AutoPtr<libMesh::TimeSolver> prev_time_solver(system->time_solver);
74 
75  libMesh::SteadySolver* steady_solver = new libMesh::SteadySolver( *(system) );
76 
77  system->time_solver = libMesh::AutoPtr<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 = prev_time_solver;
94 
95  return;
96  }
97 
99  (std::tr1::shared_ptr<libMesh::EquationSystems> equation_system,
100  MultiphysicsSystem* system,
101  const libMesh::ParameterVector & params,
102  const unsigned int time_step,
103  const libMesh::Real time )
104  {
105  for (unsigned int p=0; p != params.size(); ++p)
106  {
107  std::stringstream suffix;
108  suffix << time_step;
109 
110  std::stringstream pstr;
111  pstr << p;
112 
113  std::string filename =
114  this->_vis_output_file_prefix + "_unsteady_dRdp" +
115  pstr.str() + '.' + suffix.str();
116 
117  // Swap solution with precomputed sensitivity rhs
118  system->solution->swap(system->get_sensitivity_rhs(p));
119  equation_system->update();
120 
121  this->dump_visualization( equation_system, filename, time );
122 
123  // Now swap back and reupdate
124  system->solution->swap(system->get_sensitivity_rhs(p));
125  equation_system->update();
126  }
127  }
128 
130  ( std::tr1::shared_ptr<libMesh::EquationSystems> equation_system,
131  MultiphysicsSystem* system,
132  const unsigned int time_step,
133  const libMesh::Real time )
134  {
135  std::stringstream suffix;
136  suffix << time_step;
137 
138  const libMesh::DifferentiableQoI* raw_qoi = system->get_qoi();
139  const CompositeQoI* qoi = dynamic_cast<const CompositeQoI*>( raw_qoi );
140 
141  unsigned int n_qois = qoi->n_qois();
142 
143  for( unsigned int q = 0; q < n_qois; q++ )
144  {
145  libMesh::NumericVector<libMesh::Number>& dual_solution = system->get_adjoint_solution(q);
146 
147  const std::string& qoi_name = qoi->get_qoi(q).name();
148  std::string filename = this->_vis_output_file_prefix+"_unsteady_adjoint_"+qoi_name;
149  filename+="."+suffix.str();
150 
151  system->solution->swap( dual_solution );
152  equation_system->update();
153 
154  this->dump_visualization( equation_system, filename, time );
155 
156  // Now swap back and reupdate
157  system->solution->swap( dual_solution );
158  equation_system->update();
159  }
160  }
161 
163  (std::tr1::shared_ptr<libMesh::EquationSystems> equation_system,
164  MultiphysicsSystem* system,
165  const libMesh::ParameterVector & params,
166  const unsigned int time_step,
167  const libMesh::Real time )
168  {
169  for (unsigned int p=0; p != params.size(); ++p)
170  {
171  std::stringstream suffix;
172  suffix << time_step;
173 
174  std::stringstream pstr;
175  pstr << p;
176 
177  std::string filename =
178  this->_vis_output_file_prefix + "_unsteady_dudp" +
179  pstr.str() + '.' + suffix.str();
180 
181  // Swap solution with precomputed sensitivity solution
182  system->solution->swap(system->get_sensitivity_solution(p));
183  equation_system->update();
184 
185  this->dump_visualization( equation_system, filename, time );
186 
187  // Now swap back and reupdate
188  system->solution->swap(system->get_sensitivity_solution(p));
189  equation_system->update();
190  }
191  }
192 
193 
194 } // namespace GRINS
const std::string & name() const
Returns the name of this QoI.
Definition: qoi_base.h:143
virtual void output_adjoint(std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const unsigned int time_step, const libMesh::Real time)
const QoIBase & get_qoi(unsigned int qoi_index) const
GRINS namespace.
UnsteadyVisualization(const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
virtual void output_residual(std::tr1::shared_ptr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const unsigned int time_step, const libMesh::Real time)
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.
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)
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