GRINS-0.7.0
List of all members | Public Member Functions
GRINS::UnsteadyVisualization Class Reference

#include <unsteady_visualization.h>

Inheritance diagram for GRINS::UnsteadyVisualization:
Inheritance graph
[legend]
Collaboration diagram for GRINS::UnsteadyVisualization:
Collaboration graph
[legend]

Public Member Functions

 UnsteadyVisualization (const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
 ~UnsteadyVisualization ()
 
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)
 
virtual void output_adjoint (SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const unsigned int time_step, const libMesh::Real time)
 
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)
 
- Public Member Functions inherited from GRINS::Visualization
 Visualization (const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
 
virtual ~Visualization ()
 
void output (SharedPtr< libMesh::EquationSystems > equation_system)
 
void output (SharedPtr< libMesh::EquationSystems > equation_system, const unsigned int time_step, const libMesh::Real time)
 
void output_residual (SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system)
 
void output_residual_sensitivities (SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const libMesh::ParameterVector &params)
 
void output_adjoint (SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system)
 
void output_solution_sensitivities (SharedPtr< libMesh::EquationSystems > equation_system, GRINS::MultiphysicsSystem *system, const libMesh::ParameterVector &params)
 
void dump_visualization (SharedPtr< libMesh::EquationSystems > equation_system, const std::string &filename_prefix, const libMesh::Real time)
 

Additional Inherited Members

- Protected Attributes inherited from GRINS::Visualization
std::string _vis_output_file_prefix
 
std::vector< std::string > _output_format
 

Detailed Description

Definition at line 34 of file unsteady_visualization.h.

Constructor & Destructor Documentation

GRINS::UnsteadyVisualization::UnsteadyVisualization ( const GetPot &  input,
const libMesh::Parallel::Communicator &comm  LIBMESH_CAN_DEFAULT_TO_COMMWORLD 
)

Definition at line 42 of file unsteady_visualization.C.

44  : Visualization(input, comm)
45  {
46  return;
47  }
Visualization(const GetPot &input, const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
Definition: visualization.C:50
GRINS::UnsteadyVisualization::~UnsteadyVisualization ( )

Definition at line 49 of file unsteady_visualization.C.

50  {
51  return;
52  }

Member Function Documentation

void GRINS::UnsteadyVisualization::output_adjoint ( SharedPtr< libMesh::EquationSystems >  equation_system,
GRINS::MultiphysicsSystem system,
const unsigned int  time_step,
const libMesh::Real  time 
)
virtual

Implements GRINS::Visualization.

Definition at line 128 of file unsteady_visualization.C.

References GRINS::CompositeQoI::get_qoi(), GRINS::CompositeQoI::n_qois(), and GRINS::QoIBase::name().

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  }
void dump_visualization(SharedPtr< libMesh::EquationSystems > equation_system, const std::string &filename_prefix, const libMesh::Real time)
std::string _vis_output_file_prefix
void GRINS::UnsteadyVisualization::output_residual ( SharedPtr< libMesh::EquationSystems >  equation_system,
GRINS::MultiphysicsSystem system,
const unsigned int  time_step,
const libMesh::Real  time 
)
virtual

Implements GRINS::Visualization.

Definition at line 55 of file unsteady_visualization.C.

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  }
void dump_visualization(SharedPtr< libMesh::EquationSystems > equation_system, const std::string &filename_prefix, const libMesh::Real time)
std::string _vis_output_file_prefix
void GRINS::UnsteadyVisualization::output_residual_sensitivities ( SharedPtr< libMesh::EquationSystems >  equation_system,
GRINS::MultiphysicsSystem system,
const libMesh::ParameterVector &  params,
const unsigned int  time_step,
const libMesh::Real  time 
)
virtual

Implements GRINS::Visualization.

Definition at line 97 of file unsteady_visualization.C.

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  }
void dump_visualization(SharedPtr< libMesh::EquationSystems > equation_system, const std::string &filename_prefix, const libMesh::Real time)
std::string _vis_output_file_prefix
void GRINS::UnsteadyVisualization::output_solution_sensitivities ( SharedPtr< libMesh::EquationSystems >  equation_system,
GRINS::MultiphysicsSystem system,
const libMesh::ParameterVector &  params,
const unsigned int  time_step,
const libMesh::Real  time 
)
virtual

Implements GRINS::Visualization.

Definition at line 161 of file unsteady_visualization.C.

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  }
void dump_visualization(SharedPtr< libMesh::EquationSystems > equation_system, const std::string &filename_prefix, const libMesh::Real time)
std::string _vis_output_file_prefix

The documentation for this class was generated from the following files:

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