34 #include "libmesh/getpot.h"
35 #include "libmesh/parameter_vector.h"
36 #include "libmesh/steady_solver.h"
42 (
const GetPot& input,
43 const libMesh::Parallel::Communicator &comm )
55 ( SharedPtr<libMesh::EquationSystems> equation_system,
57 const unsigned int time_step,
58 const libMesh::Real time )
60 std::stringstream suffix;
63 std::string filename = this->_vis_output_file_prefix+
"_unsteady_residual";
65 filename+=
"."+suffix.str();
73 libMesh::TimeSolver* prev_time_solver = system->time_solver.get();
75 libMesh::SteadySolver* steady_solver =
new libMesh::SteadySolver( *(system) );
77 system->time_solver = libMesh::UniquePtr<libMesh::TimeSolver>(steady_solver);
83 system->solution->swap( *(system->rhs) );
85 equation_system->update();
87 this->dump_visualization( equation_system, filename, time );
90 system->solution->swap( *(system->rhs) );
91 equation_system->update();
93 system->time_solver.reset(prev_time_solver);
97 (SharedPtr<libMesh::EquationSystems> equation_system,
99 const libMesh::ParameterVector & params,
100 const unsigned int time_step,
101 const libMesh::Real time )
103 for (
unsigned int p=0; p != params.size(); ++p)
105 std::stringstream suffix;
108 std::stringstream pstr;
111 std::string filename =
112 this->_vis_output_file_prefix +
"_unsteady_dRdp" +
113 pstr.str() +
'.' + suffix.str();
116 system->solution->swap(system->get_sensitivity_rhs(p));
117 equation_system->update();
119 this->dump_visualization( equation_system, filename, time );
122 system->solution->swap(system->get_sensitivity_rhs(p));
123 equation_system->update();
128 ( SharedPtr<libMesh::EquationSystems> equation_system,
130 const unsigned int time_step,
131 const libMesh::Real time )
133 std::stringstream suffix;
136 const libMesh::DifferentiableQoI* raw_qoi = system->get_qoi();
139 unsigned int n_qois = qoi->
n_qois();
141 for(
unsigned int q = 0; q < n_qois; q++ )
143 libMesh::NumericVector<libMesh::Number>& dual_solution = system->get_adjoint_solution(q);
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();
149 system->solution->swap( dual_solution );
150 equation_system->update();
152 this->dump_visualization( equation_system, filename, time );
155 system->solution->swap( dual_solution );
156 equation_system->update();
161 (SharedPtr<libMesh::EquationSystems> equation_system,
163 const libMesh::ParameterVector & params,
164 const unsigned int time_step,
165 const libMesh::Real time )
167 for (
unsigned int p=0; p != params.size(); ++p)
169 std::stringstream suffix;
172 std::stringstream pstr;
175 std::string filename =
176 this->_vis_output_file_prefix +
"_unsteady_dudp" +
177 pstr.str() +
'.' + suffix.str();
180 system->solution->swap(system->get_sensitivity_solution(p));
181 equation_system->update();
183 this->dump_visualization( equation_system, filename, time );
186 system->solution->swap(system->get_sensitivity_solution(p));
187 equation_system->update();
const std::string & name() const
Returns the name of this QoI.
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 ¶ms, 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 ¶ms, const unsigned int time_step, const libMesh::Real time)
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.
virtual void assembly(bool get_residual, bool get_jacobian, bool apply_heterogeneous_constraints=false, bool apply_no_constraints=false)
Override FEMSystem::assembly.
unsigned int n_qois() const