GRINS-0.7.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, CachedValues &cache)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 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 init_variables (libMesh::FEMSystem *system)
 Initialization of variables. More...
 
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 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 side_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 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 (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_damping_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
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)
 
- 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
 
VariableIndex _fan_speed_var
 
std::string _fan_speed_var_name
 
- Protected Attributes inherited from GRINS::IncompressibleNavierStokesBase< Viscosity >
unsigned int _dim
 Physical dimension of problem. More...
 
VelocityFEVariables _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 
47  return;
48  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
template<class Mu >
GRINS::AveragedTurbine< Mu >::~AveragedTurbine ( )

Definition at line 51 of file averaged_turbine.C.

52  {
53  return;
54  }
template<class Viscosity >
GRINS::AveragedTurbine< Viscosity >::AveragedTurbine ( )
private

Member Function Documentation

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

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

Reimplemented from GRINS::Physics.

Definition at line 68 of file averaged_turbine.C.

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

59  {
60  context.get_element_fe(this->_flow_vars.u())->get_xyz();
61  context.get_element_fe(this->_flow_vars.u())->get_phi();
62 
63  return;
64  }
template<class Mu >
void GRINS::AveragedTurbine< Mu >::nonlocal_mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Mass matrix part(s) for scalar variables.

Reimplemented from GRINS::Physics.

Definition at line 289 of file averaged_turbine.C.

292  {
293  libMesh::DenseSubMatrix<libMesh::Number> &Kss =
294  context.get_elem_jacobian(this->fan_speed_var(), this->fan_speed_var()); // R_{s},{s}
295 
296  libMesh::DenseSubVector<libMesh::Number> &Fs =
297  context.get_elem_residual(this->fan_speed_var()); // R_{s}
298 
299  const libMesh::DenseSubVector<libMesh::Number> &Us =
300  context.get_elem_solution_rate(this->_fan_speed_var);
301 
302  const libMesh::Number& fan_speed = Us(0);
303 
304  Fs(0) -= this->moment_of_inertia * fan_speed;
305 
306  if (compute_jacobian)
307  {
308  Kss(0,0) -= this->moment_of_inertia * context.get_elem_solution_rate_derivative();
309  }
310 
311  return;
312  }
VariableIndex fan_speed_var() const
template<class Mu >
void GRINS::AveragedTurbine< Mu >::nonlocal_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 251 of file averaged_turbine.C.

254  {
255  libMesh::DenseSubMatrix<libMesh::Number> &Kss =
256  context.get_elem_jacobian(this->fan_speed_var(), this->fan_speed_var()); // R_{s},{s}
257 
258  libMesh::DenseSubVector<libMesh::Number> &Fs =
259  context.get_elem_residual(this->fan_speed_var()); // R_{s}
260 
261  const std::vector<libMesh::dof_id_type>& dof_indices =
262  context.get_dof_indices(this->fan_speed_var());
263 
264  const libMesh::Number fan_speed =
265  context.get_system().current_solution(dof_indices[0]);
266 
267  const libMesh::Number output_torque =
268  this->torque_function(libMesh::Point(0), fan_speed);
269 
270  Fs(0) += output_torque;
271 
272  if (compute_jacobian)
273  {
274  // FIXME: we should replace this FEM with a hook to the AD fparser stuff
275  const libMesh::Number epsilon = 1e-6;
276  const libMesh::Number output_torque_deriv =
277  (this->torque_function(libMesh::Point(0), fan_speed+epsilon) -
278  this->torque_function(libMesh::Point(0), fan_speed-epsilon)) / (2*epsilon);
279 
280  Kss(0,0) += output_torque_deriv * context.get_elem_solution_derivative();
281  }
282 
283  return;
284  }
libMesh::ParsedFunction< libMesh::Number > torque_function
VariableIndex fan_speed_var() const

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

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