GRINS-0.6.0
averaged_turbine.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2015 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 
26 // This class
27 #include "grins/averaged_turbine.h"
28 
29 // GRINS
33 
34 // libMesh
35 #include "libmesh/quadrature.h"
36 #include "libmesh/boundary_info.h"
37 
38 namespace GRINS
39 {
40 
41  template<class Mu>
42  AveragedTurbine<Mu>::AveragedTurbine( const std::string& physics_name, const GetPot& input )
43  : AveragedTurbineBase<Mu>(physics_name, input)
44  {
45  this->_ic_handler = new GenericICHandler( physics_name, input );
46 
47  return;
48  }
49 
50  template<class Mu>
52  {
53  return;
54  }
55 
56 
57  template<class Mu>
59  {
60  context.get_element_fe(this->_flow_vars.u_var())->get_xyz();
61  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
62 
63  return;
64  }
65 
66 
67  template<class Mu>
68  void AveragedTurbine<Mu>::element_time_derivative( bool compute_jacobian,
69  AssemblyContext& context,
70  CachedValues& /* cache */ )
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_var())->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_var())->get_phi();
83 
84  const std::vector<libMesh::Point>& u_qpoint =
85  context.get_element_fe(this->_flow_vars.u_var())->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_var()).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_var(), this->_flow_vars.u_var()); // R_{u},{u}
92  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.v_var()); // R_{u},{v}
93  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.u_var()); // R_{v},{u}
94  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.v_var()); // R_{v},{v}
95 
96  libMesh::DenseSubMatrix<libMesh::Number> &Kus =
97  context.get_elem_jacobian(this->_flow_vars.u_var(),
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_var()); // R_{s},{u}
102  libMesh::DenseSubMatrix<libMesh::Number> &Kvs =
103  context.get_elem_jacobian(this->_flow_vars.v_var(),
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_var()); // 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_var()); // R_{u}
122  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v_var()); // 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_var(), this->_flow_vars.w_var()); // R_{u},{w}
130  Kvw = &context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.w_var()); // R_{v},{w}
131 
132  Kwu = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.u_var()); // R_{w},{u}
133  Kwv = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.v_var()); // R_{w},{v}
134  Kww = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.w_var()); // R_{w},{w}
135  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
136 
137  Ksw = &context.get_elem_jacobian(this->fan_speed_var(), this->_flow_vars.w_var()); // R_{s},{w}
138  Kws = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->fan_speed_var()); // R_{w},{s}
139 
140  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // 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_var(), qp);
150  v = context.interior_value(this->_flow_vars.v_var(), 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_var(), 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  }
247 
248 
249 
250  template<class Mu>
252  AssemblyContext& context,
253  CachedValues& /* cache */ )
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  }
285 
286 
287 
288  template<class Mu>
289  void AveragedTurbine<Mu>::nonlocal_mass_residual( bool compute_jacobian,
290  AssemblyContext& context,
291  CachedValues& /* cache */ )
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  }
313 
314 
315 } // namespace GRINS
316 
317 // Instantiate
318 INSTANTIATE_INC_NS_SUBCLASS(AveragedTurbine);
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
virtual void nonlocal_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for scalar variables.
Base class for reading and handling initial conditions for physics classes.
GRINS namespace.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
virtual void nonlocal_mass_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Mass matrix part(s) for scalar variables.
INSTANTIATE_INC_NS_SUBCLASS(AveragedTurbine)
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.

Generated on Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1