GRINS-0.7.0
velocity_drag.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2016 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/velocity_drag.h"
28 
29 // GRINS
32 
33 // libMesh
34 #include "libmesh/quadrature.h"
35 #include "libmesh/boundary_info.h"
36 
37 namespace GRINS
38 {
39 
40  template<class Mu>
41  VelocityDrag<Mu>::VelocityDrag( const std::string& physics_name, const GetPot& input )
42  : VelocityDragBase<Mu>(physics_name, input)
43  {}
44 
45  template<class Mu>
47  {
48  context.get_element_fe(this->_flow_vars.u())->get_xyz();
49  context.get_element_fe(this->_flow_vars.u())->get_phi();
50 
51  return;
52  }
53 
54 
55  template<class Mu>
56  void VelocityDrag<Mu>::element_time_derivative( bool compute_jacobian,
57  AssemblyContext& context,
58  CachedValues& /* cache */ )
59  {
60 #ifdef GRINS_USE_GRVY_TIMERS
61  this->_timer->BeginTimer("VelocityDrag::element_time_derivative");
62 #endif
63 
64  // Element Jacobian * quadrature weights for interior integration
65  const std::vector<libMesh::Real> &JxW =
66  context.get_element_fe(this->_flow_vars.u())->get_JxW();
67 
68  // The shape functions at interior quadrature points.
69  const std::vector<std::vector<libMesh::Real> >& u_phi =
70  context.get_element_fe(this->_flow_vars.u())->get_phi();
71 
72  const std::vector<libMesh::Point>& u_qpoint =
73  context.get_element_fe(this->_flow_vars.u())->get_xyz();
74 
75  // The number of local degrees of freedom in each variable
76  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
77 
78  // The subvectors and submatrices we need to fill:
79  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
80  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // R_{u},{v}
81  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // R_{v},{u}
82  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
83 
84  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
85  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
86  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
87  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
88  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
89 
90  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
91  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
92  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
93 
94  if( this->_dim == 3 )
95  {
96  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
97  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
98 
99  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
100  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
101  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
102  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
103  }
104 
105  unsigned int n_qpoints = context.get_element_qrule().n_points();
106 
107  for (unsigned int qp=0; qp != n_qpoints; qp++)
108  {
109  // Compute the solution at the old Newton iterate.
110  libMesh::Number u, v;
111  u = context.interior_value(this->_flow_vars.u(), qp);
112  v = context.interior_value(this->_flow_vars.v(), qp);
113 
114  libMesh::NumberVectorValue U(u,v);
115  if (this->_dim == 3)
116  U(2) = context.interior_value(this->_flow_vars.w(), qp); // w
117 
118  libMesh::NumberVectorValue F;
119  libMesh::NumberTensorValue dFdU;
120  libMesh::NumberTensorValue* dFdU_ptr =
121  compute_jacobian ? &dFdU : NULL;
122  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
123  continue;
124 
125  libMesh::Real jac = JxW[qp];
126 
127  for (unsigned int i=0; i != n_u_dofs; i++)
128  {
129  const libMesh::Number jac_i = jac * u_phi[i][qp];
130 
131  Fu(i) += F(0)*jac_i;
132  Fv(i) += F(1)*jac_i;
133 
134  if( this->_dim == 3 )
135  (*Fw)(i) += F(2)*jac_i;
136 
137  if( compute_jacobian )
138  {
139  for (unsigned int j=0; j != n_u_dofs; j++)
140  {
141  const libMesh::Number jac_ij =
142  jac_i * context.get_elem_solution_derivative() *
143  u_phi[j][qp];
144  Kuu(i,j) += jac_ij * dFdU(0,0);
145  Kuv(i,j) += jac_ij * dFdU(0,1);
146  Kvu(i,j) += jac_ij * dFdU(1,0);
147  Kvv(i,j) += jac_ij * dFdU(1,1);
148 
149  if( this->_dim == 3 )
150  {
151  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
152  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
153  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
154  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
155  (*Kww)(i,j) += jac_ij * dFdU(2,2);
156  }
157  }
158  }
159  }
160  }
161 
162 
163 #ifdef GRINS_USE_GRVY_TIMERS
164  this->_timer->EndTimer("VelocityDrag::element_time_derivative");
165 #endif
166 
167  return;
168  }
169 
170 } // namespace GRINS
171 
172 // Instantiate
173 INSTANTIATE_INC_NS_SUBCLASS(VelocityDrag);
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
Definition: velocity_drag.C:46
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
Definition: velocity_drag.C:56
GRINS namespace.
INSTANTIATE_INC_NS_SUBCLASS(VelocityDrag)

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