GRINS-0.8.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-2017 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 
48  template<class Mu>
50  {
51  context.get_element_fe(this->_flow_vars.u())->get_xyz();
52  context.get_element_fe(this->_flow_vars.u())->get_phi();
53  }
54 
55 
56  template<class Mu>
58  ( bool compute_jacobian, AssemblyContext & context )
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  }
224 
225 
226 
227  template<class Mu>
229  (bool compute_jacobian, AssemblyContext & context )
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  }
261 
262 
263 
264  template<class Mu>
266  ( bool compute_jacobian, AssemblyContext & context )
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  }
288 
289 
290 } // namespace GRINS
291 
292 // Instantiate
293 INSTANTIATE_INC_NS_SUBCLASS(AveragedTurbine);
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
Base class for reading and handling initial conditions for physics classes.
GRINS namespace.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
INSTANTIATE_INC_NS_SUBCLASS(AveragedTurbine)
virtual void nonlocal_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for scalar variables.
virtual void nonlocal_mass_residual(bool compute_jacobian, AssemblyContext &context)
Mass matrix part(s) for scalar variables.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.

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