GRINS-0.8.0
velocity_penalty_adjoint_stab.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_config.h"
29 
30 // GRINS
31 #include "grins/assembly_context.h"
33 
34 // libMesh
35 #include "libmesh/getpot.h"
36 #include "libmesh/quadrature.h"
37 
38 namespace GRINS
39 {
40 
41  template<class Mu>
42  VelocityPenaltyAdjointStabilization<Mu>::VelocityPenaltyAdjointStabilization( const std::string& physics_name, const GetPot& input )
43  : VelocityPenaltyBase<Mu>(physics_name,input),
44  _stab_helper( physics_name+"StabHelper", input )
45  {
46  }
47 
48  template<class Mu>
50  {
51  return;
52  }
53 
54  template<class Mu>
56  {
57  context.get_element_fe(this->_press_var.p())->get_dphi();
58 
59  context.get_element_fe(this->_flow_vars.u())->get_xyz();
60  context.get_element_fe(this->_flow_vars.u())->get_phi();
61  context.get_element_fe(this->_flow_vars.u())->get_dphi();
62  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
63 
64  return;
65  }
66 
67  template<class Mu>
69  ( bool compute_jacobian,
70  AssemblyContext & context )
71  {
72  // The number of local degrees of freedom in each variable.
73  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
74 
75  // Element Jacobian * quadrature weights for interior integration.
76  const std::vector<libMesh::Real> &JxW =
77  context.get_element_fe(this->_flow_vars.u())->get_JxW();
78 
79  const std::vector<libMesh::Point>& u_qpoint =
80  context.get_element_fe(this->_flow_vars.u())->get_xyz();
81 
82  const std::vector<std::vector<libMesh::Real> >& u_phi =
83  context.get_element_fe(this->_flow_vars.u())->get_phi();
84 
85  const std::vector<std::vector<libMesh::RealGradient> >& u_gradphi =
86  context.get_element_fe(this->_flow_vars.u())->get_dphi();
87 
88  const std::vector<std::vector<libMesh::RealTensor> >& u_hessphi =
89  context.get_element_fe(this->_flow_vars.u())->get_d2phi();
90 
91  // Get residuals and jacobians
92  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
93  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
94  libMesh::DenseSubVector<libMesh::Number> *Fw = NULL;
95 
96  libMesh::DenseSubMatrix<libMesh::Number> &Kuu =
97  context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // J_{uu}
98  libMesh::DenseSubMatrix<libMesh::Number> &Kuv =
99  context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // J_{uv}
100  libMesh::DenseSubMatrix<libMesh::Number> &Kvu =
101  context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // J_{vu}
102  libMesh::DenseSubMatrix<libMesh::Number> &Kvv =
103  context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // J_{vv}
104 
105  libMesh::DenseSubMatrix<libMesh::Number> *Kuw = NULL;
106  libMesh::DenseSubMatrix<libMesh::Number> *Kvw = NULL;
107  libMesh::DenseSubMatrix<libMesh::Number> *Kwu = NULL;
108  libMesh::DenseSubMatrix<libMesh::Number> *Kwv = NULL;
109  libMesh::DenseSubMatrix<libMesh::Number> *Kww = NULL;
110 
111  if(this->_flow_vars.dim() == 3)
112  {
113  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
114  Kuw = &context.get_elem_jacobian
115  (this->_flow_vars.u(), this->_flow_vars.w()); // J_{uw}
116  Kvw = &context.get_elem_jacobian
117  (this->_flow_vars.v(), this->_flow_vars.w()); // J_{vw}
118  Kwu = &context.get_elem_jacobian
119  (this->_flow_vars.w(), this->_flow_vars.u()); // J_{wu}
120  Kwv = &context.get_elem_jacobian
121  (this->_flow_vars.w(), this->_flow_vars.v()); // J_{wv}
122  Kww = &context.get_elem_jacobian
123  (this->_flow_vars.w(), this->_flow_vars.w()); // J_{ww}
124  }
125 
126  // Now we will build the element Jacobian and residual.
127  // Constructing the residual requires the solution and its
128  // gradient from the previous timestep. This must be
129  // calculated at each quadrature point by summing the
130  // solution degree-of-freedom values by the appropriate
131  // weight functions.
132  unsigned int n_qpoints = context.get_element_qrule().n_points();
133 
134  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
135 
136  for (unsigned int qp=0; qp != n_qpoints; qp++)
137  {
138  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
139  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
140 
141  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
142  context.interior_value( this->_flow_vars.v(), qp ) );
143  if( this->_flow_vars.dim() == 3 )
144  {
145  U(2) = context.interior_value( this->_flow_vars.w(), qp );
146  }
147 
148  // Compute the viscosity at this qp
149  libMesh::Real mu_qp = this->_mu(context, qp);
150 
151  libMesh::Real tau_M;
152  libMesh::Real d_tau_M_d_rho;
153  libMesh::Gradient d_tau_M_dU;
154 
155  if (compute_jacobian)
156  this->_stab_helper.compute_tau_momentum_and_derivs
157  ( context, qp, g, G, this->_rho, U, mu_qp,
158  tau_M, d_tau_M_d_rho, d_tau_M_dU,
159  this->_is_steady );
160  else
161  tau_M = this->_stab_helper.compute_tau_momentum
162  ( context, qp, g, G, this->_rho, U, mu_qp,
163  this->_is_steady );
164 
165  libMesh::NumberVectorValue F;
166  libMesh::NumberTensorValue dFdU;
167  libMesh::NumberTensorValue* dFdU_ptr =
168  compute_jacobian ? &dFdU : NULL;
169  if (!this->compute_force(u_qpoint[qp], context, U, F, dFdU_ptr))
170  continue;
171 
172  for (unsigned int i=0; i != n_u_dofs; i++)
173  {
174  libMesh::Real test_func = this->_rho*U*u_gradphi[i][qp] +
175  mu_qp*( u_hessphi[i][qp](0,0) + u_hessphi[i][qp](1,1) + u_hessphi[i][qp](2,2) );
176  Fu(i) += tau_M*F(0)*test_func*JxW[qp];
177 
178  Fv(i) += tau_M*F(1)*test_func*JxW[qp];
179 
180  if (this->_flow_vars.dim() == 3)
181  {
182  (*Fw)(i) += tau_M*F(2)*test_func*JxW[qp];
183  }
184 
185  if (compute_jacobian)
186  {
187  libMesh::Gradient d_test_func_dU = this->_rho*u_gradphi[i][qp];
188 
189  for (unsigned int j=0; j != n_u_dofs; ++j)
190  {
191  Kuu(i,j) += tau_M*F(0)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
192  Kuu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
193  Kuu(i,j) += tau_M*dFdU(0,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
194  Kuv(i,j) += tau_M*F(0)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
195  Kuv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
196  Kuv(i,j) += tau_M*dFdU(0,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
197  Kvu(i,j) += tau_M*F(1)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
198  Kvu(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
199  Kvu(i,j) += tau_M*dFdU(1,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
200  Kvv(i,j) += tau_M*F(1)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
201  Kvv(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
202  Kvv(i,j) += tau_M*dFdU(1,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
203  }
204 
205  if (this->_flow_vars.dim() == 3)
206  {
207  for (unsigned int j=0; j != n_u_dofs; ++j)
208  {
209  (*Kuw)(i,j) += tau_M*F(0)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
210  (*Kuw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(0)*test_func*JxW[qp]*context.get_elem_solution_derivative();
211  (*Kuw)(i,j) += tau_M*dFdU(0,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
212  (*Kvw)(i,j) += tau_M*F(1)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
213  (*Kvw)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(1)*test_func*JxW[qp]*context.get_elem_solution_derivative();
214  (*Kvw)(i,j) += tau_M*dFdU(1,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
215  (*Kwu)(i,j) += tau_M*F(2)*d_test_func_dU(0)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
216  (*Kwu)(i,j) += d_tau_M_dU(0)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
217  (*Kwu)(i,j) += tau_M*dFdU(2,0)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
218  (*Kwv)(i,j) += tau_M*F(2)*d_test_func_dU(1)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
219  (*Kwv)(i,j) += d_tau_M_dU(1)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
220  (*Kwv)(i,j) += tau_M*dFdU(2,1)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
221  (*Kww)(i,j) += tau_M*F(2)*d_test_func_dU(2)*u_phi[j][qp]*JxW[qp]*context.get_elem_solution_derivative();
222  (*Kww)(i,j) += d_tau_M_dU(2)*u_phi[j][qp]*F(2)*test_func*JxW[qp]*context.get_elem_solution_derivative();
223  (*Kww)(i,j) += tau_M*dFdU(2,2)*u_phi[j][qp]*test_func*JxW[qp]*context.get_elem_solution_derivative();
224  }
225  }
226 
227  } // End compute_jacobian check
228 
229  } // End i dof loop
230  } // End quadrature loop
231  }
232 
233  template<class Mu>
235  ( bool compute_jacobian, AssemblyContext & context )
236  {
237  // The number of local degrees of freedom in each variable.
238  const unsigned int n_p_dofs = context.get_dof_indices(this->_press_var.p()).size();
239  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
240 
241  // Element Jacobian * quadrature weights for interior integration.
242  const std::vector<libMesh::Real> &JxW =
243  context.get_element_fe(this->_flow_vars.u())->get_JxW();
244 
245  const std::vector<libMesh::Point>& u_qpoint =
246  context.get_element_fe(this->_flow_vars.u())->get_xyz();
247 
248  const std::vector<std::vector<libMesh::Real> >& u_phi =
249  context.get_element_fe(this->_flow_vars.u())->get_phi();
250 
251  const std::vector<std::vector<libMesh::RealGradient> >& p_dphi =
252  context.get_element_fe(this->_press_var.p())->get_dphi();
253 
254  libMesh::DenseSubVector<libMesh::Number> &Fp = context.get_elem_residual(this->_press_var.p()); // R_{p}
255 
256  libMesh::DenseSubMatrix<libMesh::Number> &Kpu =
257  context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.u()); // J_{pu}
258  libMesh::DenseSubMatrix<libMesh::Number> &Kpv =
259  context.get_elem_jacobian(this->_press_var.p(), this->_flow_vars.v()); // J_{pv}
260  libMesh::DenseSubMatrix<libMesh::Number> *Kpw = NULL;
261 
262  if(this->_flow_vars.dim() == 3)
263  {
264  Kpw = &context.get_elem_jacobian
265  (this->_press_var.p(), this->_flow_vars.w()); // J_{pw}
266  }
267 
268  // Now we will build the element Jacobian and residual.
269  // Constructing the residual requires the solution and its
270  // gradient from the previous timestep. This must be
271  // calculated at each quadrature point by summing the
272  // solution degree-of-freedom values by the appropriate
273  // weight functions.
274  unsigned int n_qpoints = context.get_element_qrule().n_points();
275 
276  libMesh::FEBase* fe = context.get_element_fe(this->_flow_vars.u());
277 
278  for (unsigned int qp=0; qp != n_qpoints; qp++)
279  {
280  libMesh::RealGradient g = this->_stab_helper.compute_g( fe, context, qp );
281  libMesh::RealTensor G = this->_stab_helper.compute_G( fe, context, qp );
282 
283  libMesh::RealGradient U( context.interior_value( this->_flow_vars.u(), qp ),
284  context.interior_value( this->_flow_vars.v(), qp ) );
285  if( this->_flow_vars.dim() == 3 )
286  {
287  U(2) = context.interior_value( this->_flow_vars.w(), qp );
288  }
289 
290  // Compute the viscosity at this qp
291  libMesh::Real mu_qp = this->_mu(context, qp);
292 
293  libMesh::Real tau_M;
294  libMesh::Real d_tau_M_d_rho;
295  libMesh::Gradient d_tau_M_dU;
296 
297  if (compute_jacobian)
298  this->_stab_helper.compute_tau_momentum_and_derivs
299  ( context, qp, g, G, this->_rho, U, mu_qp,
300  tau_M, d_tau_M_d_rho, d_tau_M_dU,
301  this->_is_steady );
302  else
303  tau_M = this->_stab_helper.compute_tau_momentum
304  ( context, qp, g, G, this->_rho, U, mu_qp,
305  this->_is_steady );
306 
307  libMesh::NumberVectorValue F;
308  libMesh::NumberTensorValue dFdU;
309  libMesh::NumberTensorValue* dFdU_ptr =
310  compute_jacobian ? &dFdU : NULL;
311  if (!this->compute_force(u_qpoint[qp], context, U, F, dFdU_ptr))
312  continue;
313 
314  // First, an i-loop over the velocity degrees of freedom.
315  // We know that n_u_dofs == n_v_dofs so we can compute contributions
316  // for both at the same time.
317  for (unsigned int i=0; i != n_p_dofs; i++)
318  {
319  Fp(i) += -tau_M*F*p_dphi[i][qp]*JxW[qp];
320 
321  if (compute_jacobian)
322  {
323  for (unsigned int j=0; j != n_u_dofs; ++j)
324  {
325  Kpu(i,j) += -d_tau_M_dU(0)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
326  Kpv(i,j) += -d_tau_M_dU(1)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
327  for (unsigned int d=0; d != 3; ++d)
328  {
329  Kpu(i,j) += -tau_M*dFdU(d,0)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
330  Kpv(i,j) += -tau_M*dFdU(d,1)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
331  }
332  }
333  if( this->_flow_vars.dim() == 3 )
334  for (unsigned int j=0; j != n_u_dofs; ++j)
335  {
336  (*Kpw)(i,j) += -d_tau_M_dU(2)*u_phi[j][qp]*F*p_dphi[i][qp]*JxW[qp]*context.get_elem_solution_derivative();
337  for (unsigned int d=0; d != 3; ++d)
338  {
339  (*Kpw)(i,j) += -tau_M*dFdU(d,2)*u_phi[j][qp]*p_dphi[i][qp](d)*JxW[qp]*context.get_elem_solution_derivative();
340  }
341  }
342  }
343  }
344  } // End quadrature loop
345  }
346 
347 } // namespace GRINS
348 
349 // Instantiate
350 INSTANTIATE_INC_NS_SUBCLASS(VelocityPenaltyAdjointStabilization);
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
GRINS namespace.
INSTANTIATE_INC_NS_SUBCLASS(VelocityPenaltyAdjointStabilization)
virtual void element_constraint(bool compute_jacobian, AssemblyContext &context)
Constraint part(s) of physics for element interiors.

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