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

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

#include <averaged_turbine.h>

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

Public Member Functions

 AveragedTurbine (const std::string &physics_name, const GetPot &input)
 
 ~AveragedTurbine ()
 
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 nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context)
 Mass matrix part(s) for scalar variables. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for scalar variables. More...
 
- Public Member Functions inherited from GRINS::AveragedTurbineBase< Viscosity >
 AveragedTurbineBase (const std::string &physics_name, const GetPot &input)
 
 ~AveragedTurbineBase ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets turbine_speed and velocity variables to be time-evolving. More...
 
bool compute_force (const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::Number s, libMesh::NumberVectorValue &U_B_1, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL, libMesh::NumberVectorValue *dFds=NULL)
 
VariableIndex fan_speed_var () const
 
- 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 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 register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. 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 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...
 
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 &)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
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

 AveragedTurbine ()
 

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::AveragedTurbineBase< 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 > torque_function
 
libMesh::Number moment_of_inertia
 
libMesh::Number initial_speed
 
libMesh::ParsedFunction< libMesh::Number > chord_function
 
libMesh::ParsedFunction< libMesh::Number > area_swept_function
 
libMesh::ParsedFunction< libMesh::Number > aoa_function
 
ScalarVariable_var
 
- 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::AveragedTurbine< Viscosity >

Physics class for spatially-averaged turbine.

Definition at line 52 of file averaged_turbine.h.

Constructor & Destructor Documentation

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

Definition at line 42 of file averaged_turbine.C.

References GRINS::Physics::_ic_handler.

43  : AveragedTurbineBase<Mu>(physics_name, input)
44  {
45  this->_ic_handler = new GenericICHandler( physics_name, input );
46  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
template<class Viscosity >
GRINS::AveragedTurbine< Viscosity >::~AveragedTurbine ( )
inline

Definition at line 58 of file averaged_turbine.h.

58 {};
template<class Viscosity >
GRINS::AveragedTurbine< Viscosity >::AveragedTurbine ( )
private

Member Function Documentation

template<class Mu >
void GRINS::AveragedTurbine< Mu >::element_time_derivative ( bool  ,
AssemblyContext  
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 58 of file averaged_turbine.C.

59  {
60  // Element Jacobian * quadrature weights for interior integration
61  const std::vector<libMesh::Real> &JxW =
62  context.get_element_fe(this->_flow_vars.u())->get_JxW();
63 
64  // The shape functions at interior quadrature points.
65  const std::vector<std::vector<libMesh::Real> >& u_phi =
66  context.get_element_fe(this->_flow_vars.u())->get_phi();
67 
68  const std::vector<libMesh::Point>& u_qpoint =
69  context.get_element_fe(this->_flow_vars.u())->get_xyz();
70 
71  // The number of local degrees of freedom in each variable
72  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
73 
74  // The subvectors and submatrices we need to fill:
75  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
76  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // R_{u},{v}
77  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // R_{v},{u}
78  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
79 
80  libMesh::DenseSubMatrix<libMesh::Number> &Kus =
81  context.get_elem_jacobian(this->_flow_vars.u(),
82  this->fan_speed_var()); // R_{u},{s}
83  libMesh::DenseSubMatrix<libMesh::Number> &Ksu =
84  context.get_elem_jacobian(this->fan_speed_var(),
85  this->_flow_vars.u()); // R_{s},{u}
86  libMesh::DenseSubMatrix<libMesh::Number> &Kvs =
87  context.get_elem_jacobian(this->_flow_vars.v(),
88  this->fan_speed_var()); // R_{v},{s}
89  libMesh::DenseSubMatrix<libMesh::Number> &Ksv =
90  context.get_elem_jacobian(this->fan_speed_var(),
91  this->_flow_vars.v()); // R_{s},{v}
92  libMesh::DenseSubMatrix<libMesh::Number> &Kss =
93  context.get_elem_jacobian(this->fan_speed_var(),
94  this->fan_speed_var()); // R_{s},{s}
95 
96  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
97  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
98  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
99  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
100  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
101 
102  libMesh::DenseSubMatrix<libMesh::Number>* Ksw = NULL;
103  libMesh::DenseSubMatrix<libMesh::Number>* Kws = NULL;
104 
105  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
106  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
107  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
108 
109  libMesh::DenseSubVector<libMesh::Number> &Fs = context.get_elem_residual(this->fan_speed_var()); // R_{s}
110 
111  if( this->_flow_vars.dim() == 3 )
112  {
113  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
114  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
115 
116  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
117  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
118  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
119  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
120 
121  Ksw = &context.get_elem_jacobian(this->fan_speed_var(), this->_flow_vars.w()); // R_{s},{w}
122  Kws = &context.get_elem_jacobian(this->_flow_vars.w(), this->fan_speed_var()); // R_{w},{s}
123 
124  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
125  }
126 
127  unsigned int n_qpoints = context.get_element_qrule().n_points();
128 
129  for (unsigned int qp=0; qp != n_qpoints; qp++)
130  {
131  // Compute the solution at the old Newton iterate.
132  libMesh::Number u, v, s;
133  u = context.interior_value(this->_flow_vars.u(), qp);
134  v = context.interior_value(this->_flow_vars.v(), qp);
135  s = context.interior_value(this->fan_speed_var(), qp);
136 
137  libMesh::NumberVectorValue U(u,v);
138  if (this->_flow_vars.dim() == 3)
139  U(2) = context.interior_value(this->_flow_vars.w(), qp); // w
140 
141  libMesh::NumberVectorValue U_B_1;
142  libMesh::NumberVectorValue F;
143  libMesh::NumberTensorValue dFdU;
144  libMesh::NumberTensorValue* dFdU_ptr =
145  compute_jacobian ? &dFdU : NULL;
146  libMesh::NumberVectorValue dFds;
147  libMesh::NumberVectorValue* dFds_ptr =
148  compute_jacobian ? &dFds : NULL;
149  if (!this->compute_force(u_qpoint[qp], context.time, U, s,
150  U_B_1, F, dFdU_ptr, dFds_ptr))
151  continue;
152 
153  libMesh::Real jac = JxW[qp];
154 
155  // Using this dot product to derive torque *depends* on s=1
156  // and U_B_1 corresponding to 1 rad/sec base velocity; this
157  // means that the length of U_B_1 is equal to radius.
158 
159  // F is the force on the air, so *negative* F is the force on
160  // the turbine.
161  Fs(0) -= U_B_1 * F * jac;
162 
163  if (compute_jacobian)
164  {
165  Kss(0,0) -= U_B_1 * dFds * jac;
166 
167  for (unsigned int j=0; j != n_u_dofs; j++)
168  {
169  libMesh::Real jac_j = JxW[qp] * u_phi[j][qp];
170 
171  for (unsigned int d=0; d != 3; ++d)
172  {
173  Ksu(0,j) -= jac_j * U_B_1(d) * dFdU(d,0);
174  Ksv(0,j) -= jac_j * U_B_1(d) * dFdU(d,1);
175  }
176 
177  if (this->_flow_vars.dim() == 3)
178  {
179  for (unsigned int d=0; d != 3; ++d)
180  (*Ksw)(0,j) -= jac_j * U_B_1(d) * dFdU(d,2);
181  }
182 
183  } // End j dof loop
184  }
185 
186  for (unsigned int i=0; i != n_u_dofs; i++)
187  {
188  const libMesh::Number jac_i = jac * u_phi[i][qp];
189 
190  Fu(i) += F(0)*jac_i;
191  Fv(i) += F(1)*jac_i;
192 
193  if( this->_flow_vars.dim() == 3 )
194  (*Fw)(i) += F(2)*jac_i;
195 
196  if( compute_jacobian )
197  {
198  Kus(i,0) += dFds(0) * jac_i;
199  Kvs(i,0) += dFds(1) * jac_i;
200  if( this->_flow_vars.dim() == 3 )
201  (*Kws)(i,0) += dFds(2) * jac_i;
202 
203  for (unsigned int j=0; j != n_u_dofs; j++)
204  {
205  const libMesh::Number jac_ij = jac_i * u_phi[j][qp];
206  Kuu(i,j) += jac_ij * dFdU(0,0);
207  Kuv(i,j) += jac_ij * dFdU(0,1);
208  Kvu(i,j) += jac_ij * dFdU(1,0);
209  Kvv(i,j) += jac_ij * dFdU(1,1);
210 
211  if( this->_flow_vars.dim() == 3 )
212  {
213  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
214  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
215  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
216  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
217  (*Kww)(i,j) += jac_ij * dFdU(2,2);
218  }
219  }
220  }
221  }
222  }
223  }
unsigned int dim() const
Number of components.
bool compute_force(const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::Number s, libMesh::NumberVectorValue &U_B_1, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL, libMesh::NumberVectorValue *dFds=NULL)
VariableIndex fan_speed_var() const
template<class Mu >
void GRINS::AveragedTurbine< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 49 of file averaged_turbine.C.

50  {
51  context.get_element_fe(this->_flow_vars.u())->get_xyz();
52  context.get_element_fe(this->_flow_vars.u())->get_phi();
53  }
template<class Mu >
void GRINS::AveragedTurbine< Mu >::nonlocal_mass_residual ( bool  ,
AssemblyContext  
)
virtual

Mass matrix part(s) for scalar variables.

Reimplemented from GRINS::Physics.

Definition at line 266 of file averaged_turbine.C.

267  {
268  libMesh::DenseSubMatrix<libMesh::Number> &Kss =
269  context.get_elem_jacobian(this->fan_speed_var(), this->fan_speed_var()); // R_{s},{s}
270 
271  libMesh::DenseSubVector<libMesh::Number> &Fs =
272  context.get_elem_residual(this->fan_speed_var()); // R_{s}
273 
274  const libMesh::DenseSubVector<libMesh::Number> &Us =
275  context.get_elem_solution_rate(this->fan_speed_var());
276 
277  const libMesh::Number& fan_speed = Us(0);
278 
279  Fs(0) -= this->moment_of_inertia * fan_speed;
280 
281  if (compute_jacobian)
282  {
283  Kss(0,0) -= this->moment_of_inertia * context.get_elem_solution_rate_derivative();
284  }
285 
286  return;
287  }
VariableIndex fan_speed_var() const
template<class Mu >
void GRINS::AveragedTurbine< Mu >::nonlocal_time_derivative ( bool  ,
AssemblyContext  
)
virtual

Time dependent part(s) of physics for scalar variables.

Reimplemented from GRINS::Physics.

Definition at line 229 of file averaged_turbine.C.

230  {
231  libMesh::DenseSubMatrix<libMesh::Number> &Kss =
232  context.get_elem_jacobian(this->fan_speed_var(), this->fan_speed_var()); // R_{s},{s}
233 
234  libMesh::DenseSubVector<libMesh::Number> &Fs =
235  context.get_elem_residual(this->fan_speed_var()); // R_{s}
236 
237  const std::vector<libMesh::dof_id_type>& dof_indices =
238  context.get_dof_indices(this->fan_speed_var());
239 
240  const libMesh::Number fan_speed =
241  context.get_system().current_solution(dof_indices[0]);
242 
243  const libMesh::Number output_torque =
244  this->torque_function(libMesh::Point(0), fan_speed);
245 
246  Fs(0) += output_torque;
247 
248  if (compute_jacobian)
249  {
250  // FIXME: we should replace this FEM with a hook to the AD fparser stuff
251  const libMesh::Number epsilon = 1e-6;
252  const libMesh::Number output_torque_deriv =
253  (this->torque_function(libMesh::Point(0), fan_speed+epsilon) -
254  this->torque_function(libMesh::Point(0), fan_speed-epsilon)) / (2*epsilon);
255 
256  Kss(0,0) += output_torque_deriv * context.get_elem_solution_derivative();
257  }
258 
259  return;
260  }
libMesh::ParsedFunction< libMesh::Number > torque_function
VariableIndex fan_speed_var() const

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