GRINS-0.8.0
List of all members | Public Member Functions | Private Member Functions | Private Attributes
GRINS::AveragedFan< Viscosity > Class Template Reference

Physics class for spatially-averaged fan. More...

#include <averaged_fan.h>

Inheritance diagram for GRINS::AveragedFan< Viscosity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::AveragedFan< Viscosity >:
Collaboration graph
[legend]

Public Member Functions

 AveragedFan (const std::string &physics_name, const GetPot &input)
 
 ~AveragedFan ()
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register postprocessing variables for visualization output. More...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 Compute value of postprocessed quantities at libMesh::Point. More...
 
- Public Member Functions inherited from GRINS::AveragedFanBase< Viscosity >
 AveragedFanBase (const std::string &physics_name, const GetPot &input)
 
 ~AveragedFanBase ()
 
bool compute_force (const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
 
- Public Member Functions inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
 IncompressibleNavierStokesBase (const std::string &my_physics_name, const std::string &core_physics_name, const GetPot &input)
 
 ~IncompressibleNavierStokesBase ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
libMesh::Real get_viscosity_value (AssemblyContext &context, unsigned int qp) const
 
- Public Member Functions inherited from GRINS::Physics
 Physics (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~Physics ()
 
virtual void init_variables (libMesh::FEMSystem *)
 Initialize variables for this physics. More...
 
virtual bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Private Member Functions

 AveragedFan ()
 

Private Attributes

unsigned int _base_velocity_x_index
 Index from registering this postprocessed quantity. More...
 
unsigned int _base_velocity_y_index
 Index from registering this postprocessed quantity. More...
 
unsigned int _base_velocity_z_index
 Index from registering this postprocessed quantity. More...
 

Additional Inherited Members

- Static Public Member Functions inherited from GRINS::Physics
static void set_is_axisymmetric (bool is_axisymmetric)
 Set whether we should treat the problem as axisymmetric. More...
 
static bool is_axisymmetric ()
 
- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 
- Protected Member Functions inherited from GRINS::Physics
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe (const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
 
void parse_enabled_subdomains (const GetPot &input, const std::string &physics_name)
 
void check_var_subdomain_consistency (const FEVariablesBase &var) const
 Check that var is enabled on at least the subdomains this Physics is. More...
 
- Protected Attributes inherited from GRINS::AveragedFanBase< Viscosity >
libMesh::ParsedFunction< libMesh::Number > base_velocity_function
 
libMesh::ParsedFunction< libMesh::Number > local_vertical_function
 
libMesh::ParsedFunction< libMesh::Number > lift_function
 
libMesh::ParsedFunction< libMesh::Number > drag_function
 
libMesh::ParsedFunction< libMesh::Number > chord_function
 
libMesh::ParsedFunction< libMesh::Number > area_swept_function
 
libMesh::ParsedFunction< libMesh::Number > aoa_function
 
- Protected Attributes inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
VelocityVariable_flow_vars
 
PressureFEVariable_press_var
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
Viscosity _mu
 Viscosity object. More...
 
- Protected Attributes inherited from GRINS::Physics
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 
- Static Protected Attributes inherited from GRINS::Physics
static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 
static bool _is_axisymmetric = false
 Caches whether we are solving an axisymmetric problem or not. More...
 

Detailed Description

template<class Viscosity>
class GRINS::AveragedFan< Viscosity >

Physics class for spatially-averaged fan.

Definition at line 50 of file averaged_fan.h.

Constructor & Destructor Documentation

template<class Mu >
GRINS::AveragedFan< Mu >::AveragedFan ( const std::string &  physics_name,
const GetPot &  input 
)

Definition at line 43 of file averaged_fan.C.

44  : AveragedFanBase<Mu>(physics_name, input),
48  {
49  return;
50  }
unsigned int _base_velocity_y_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:85
unsigned int _base_velocity_x_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:82
unsigned int _base_velocity_z_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:88
template<class Mu >
GRINS::AveragedFan< Mu >::~AveragedFan ( )

Definition at line 53 of file averaged_fan.C.

54  {
55  return;
56  }
template<class Viscosity >
GRINS::AveragedFan< Viscosity >::AveragedFan ( )
private

Member Function Documentation

template<class Mu >
void GRINS::AveragedFan< Mu >::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtual

Compute value of postprocessed quantities at libMesh::Point.

Reimplemented from GRINS::Physics.

Definition at line 204 of file averaged_fan.C.

208  {
209  libMesh::DenseVector<libMesh::Number> output_vec(3);
210 
211  if( quantity_index == this->_base_velocity_x_index )
212  {
213  this->base_velocity_function(point, context.time, output_vec);
214  value = output_vec(0);
215  }
216  else if( quantity_index == this->_base_velocity_y_index )
217  {
218  this->base_velocity_function(point, context.time, output_vec);
219  value = output_vec(1);
220  }
221  else if( quantity_index == this->_base_velocity_z_index )
222  {
223  this->base_velocity_function(point, context.time, output_vec);
224  value = output_vec(2);
225  }
226  }
unsigned int _base_velocity_y_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:85
unsigned int _base_velocity_x_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:82
libMesh::ParsedFunction< libMesh::Number > base_velocity_function
unsigned int _base_velocity_z_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:88
template<class Mu >
void GRINS::AveragedFan< Mu >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

Time dependent part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 101 of file averaged_fan.C.

103  {
104  // Element Jacobian * quadrature weights for interior integration
105  const std::vector<libMesh::Real> &JxW =
106  context.get_element_fe(this->_flow_vars.u())->get_JxW();
107 
108  // The shape functions at interior quadrature points.
109  const std::vector<std::vector<libMesh::Real> >& u_phi =
110  context.get_element_fe(this->_flow_vars.u())->get_phi();
111 
112  const std::vector<libMesh::Point>& u_qpoint =
113  context.get_element_fe(this->_flow_vars.u())->get_xyz();
114 
115  // The number of local degrees of freedom in each variable
116  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
117 
118  // The subvectors and submatrices we need to fill:
119  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
120  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // R_{u},{v}
121  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // R_{v},{u}
122  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
123 
124  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
125  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
126  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
127  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
128  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
129 
130  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
131  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
132  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
133 
134  if( this->_flow_vars.dim() == 3 )
135  {
136  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
137  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
138 
139  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
140  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
141  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
142  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
143  }
144 
145  unsigned int n_qpoints = context.get_element_qrule().n_points();
146 
147  for (unsigned int qp=0; qp != n_qpoints; qp++)
148  {
149  // Compute the solution at the old Newton iterate.
150  libMesh::Number u, v;
151  u = context.interior_value(this->_flow_vars.u(), qp);
152  v = context.interior_value(this->_flow_vars.v(), qp);
153 
154  libMesh::NumberVectorValue U(u,v);
155  if (this->_flow_vars.dim() == 3)
156  U(2) = context.interior_value(this->_flow_vars.w(), qp); // w
157 
158  libMesh::NumberVectorValue F;
159  libMesh::NumberTensorValue dFdU;
160  libMesh::NumberTensorValue* dFdU_ptr =
161  compute_jacobian ? &dFdU : NULL;
162  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
163  continue;
164 
165  libMesh::Real jac = JxW[qp];
166 
167  for (unsigned int i=0; i != n_u_dofs; i++)
168  {
169  const libMesh::Number jac_i = jac * u_phi[i][qp];
170 
171  Fu(i) += F(0)*jac_i;
172  Fv(i) += F(1)*jac_i;
173 
174  if( this->_flow_vars.dim() == 3 )
175  (*Fw)(i) += F(2)*jac_i;
176 
177  if( compute_jacobian )
178  {
179  for (unsigned int j=0; j != n_u_dofs; j++)
180  {
181  const libMesh::Number jac_ij =
182  jac_i * context.get_elem_solution_derivative() *
183  u_phi[j][qp];
184  Kuu(i,j) += jac_ij * dFdU(0,0);
185  Kuv(i,j) += jac_ij * dFdU(0,1);
186  Kvu(i,j) += jac_ij * dFdU(1,0);
187  Kvv(i,j) += jac_ij * dFdU(1,1);
188 
189  if( this->_flow_vars.dim() == 3 )
190  {
191  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
192  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
193  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
194  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
195  (*Kww)(i,j) += jac_ij * dFdU(2,2);
196  }
197  }
198  }
199  }
200  }
201  }
bool compute_force(const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
unsigned int dim() const
Number of components.
template<class Mu >
void GRINS::AveragedFan< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 60 of file averaged_fan.C.

61  {
62  context.get_element_fe(this->_flow_vars.u())->get_xyz();
63  context.get_element_fe(this->_flow_vars.u())->get_phi();
64 
65  return;
66  }
template<class Mu >
void GRINS::AveragedFan< Mu >::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtual

Register postprocessing variables for visualization output.

Reimplemented from GRINS::Physics.

Definition at line 69 of file averaged_fan.C.

References GRINS::PostProcessedQuantities< NumericType >::register_quantity().

71  {
72  std::string section = "Physics/"+this->_physics_name+"/output_vars";
73 
74  if( input.have_variable(section) )
75  {
76  unsigned int n_vars = input.vector_variable_size(section);
77 
78  for( unsigned int v = 0; v < n_vars; v++ )
79  {
80  std::string name = input(section,"DIE!",v);
81 
82  if( name == std::string("base_velocity") )
83  {
84  this->_base_velocity_x_index = postprocessing.register_quantity("base_vel_x");
85  this->_base_velocity_y_index = postprocessing.register_quantity("base_vel_y");
86  this->_base_velocity_z_index = postprocessing.register_quantity("base_vel_z" );
87  }
88  else
89  {
90  std::cerr << "Error: Invalid output_vars value for "+this->_physics_name << std::endl
91  << " Found " << name << std::endl
92  << " Acceptable values are: base_velocity" << std::endl;
93  libmesh_error();
94  }
95  }
96  }
97  }
unsigned int _base_velocity_y_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:85
unsigned int _base_velocity_x_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:82
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:256
unsigned int _base_velocity_z_index
Index from registering this postprocessed quantity.
Definition: averaged_fan.h:88

Member Data Documentation

template<class Viscosity >
unsigned int GRINS::AveragedFan< Viscosity >::_base_velocity_x_index
private

Index from registering this postprocessed quantity.

Definition at line 82 of file averaged_fan.h.

template<class Viscosity >
unsigned int GRINS::AveragedFan< Viscosity >::_base_velocity_y_index
private

Index from registering this postprocessed quantity.

Definition at line 85 of file averaged_fan.h.

template<class Viscosity >
unsigned int GRINS::AveragedFan< Viscosity >::_base_velocity_z_index
private

Index from registering this postprocessed quantity.

Definition at line 88 of file averaged_fan.h.


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

Generated on Tue Dec 19 2017 12:47:30 for GRINS-0.8.0 by  doxygen 1.8.9.1