GRINS-0.8.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-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/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 
52 
53  template<class Mu>
55  ( bool compute_jacobian,
56  AssemblyContext & context )
57  {
58  // Element Jacobian * quadrature weights for interior integration
59  const std::vector<libMesh::Real> &JxW =
60  context.get_element_fe(this->_flow_vars.u())->get_JxW();
61 
62  // The shape functions at interior quadrature points.
63  const std::vector<std::vector<libMesh::Real> >& u_phi =
64  context.get_element_fe(this->_flow_vars.u())->get_phi();
65 
66  const std::vector<libMesh::Point>& u_qpoint =
67  context.get_element_fe(this->_flow_vars.u())->get_xyz();
68 
69  // The number of local degrees of freedom in each variable
70  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
71 
72  // The subvectors and submatrices we need to fill:
73  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
74  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // R_{u},{v}
75  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // R_{v},{u}
76  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
77 
78  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
79  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
80  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
81  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
82  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
83 
84  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
85  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
86  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
87 
88  if( this->_flow_vars.dim() == 3 )
89  {
90  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
91  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
92 
93  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
94  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
95  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
96  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
97  }
98 
99  unsigned int n_qpoints = context.get_element_qrule().n_points();
100 
101  for (unsigned int qp=0; qp != n_qpoints; qp++)
102  {
103  // Compute the solution at the old Newton iterate.
104  libMesh::Number u, v;
105  u = context.interior_value(this->_flow_vars.u(), qp);
106  v = context.interior_value(this->_flow_vars.v(), qp);
107 
108  libMesh::NumberVectorValue U(u,v);
109  if (this->_flow_vars.dim() == 3)
110  U(2) = context.interior_value(this->_flow_vars.w(), qp); // w
111 
112  libMesh::NumberVectorValue F;
113  libMesh::NumberTensorValue dFdU;
114  libMesh::NumberTensorValue* dFdU_ptr =
115  compute_jacobian ? &dFdU : NULL;
116  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
117  continue;
118 
119  libMesh::Real jac = JxW[qp];
120 
121  for (unsigned int i=0; i != n_u_dofs; i++)
122  {
123  const libMesh::Number jac_i = jac * u_phi[i][qp];
124 
125  Fu(i) += F(0)*jac_i;
126  Fv(i) += F(1)*jac_i;
127 
128  if( this->_flow_vars.dim() == 3 )
129  (*Fw)(i) += F(2)*jac_i;
130 
131  if( compute_jacobian )
132  {
133  for (unsigned int j=0; j != n_u_dofs; j++)
134  {
135  const libMesh::Number jac_ij =
136  jac_i * context.get_elem_solution_derivative() *
137  u_phi[j][qp];
138  Kuu(i,j) += jac_ij * dFdU(0,0);
139  Kuv(i,j) += jac_ij * dFdU(0,1);
140  Kvu(i,j) += jac_ij * dFdU(1,0);
141  Kvv(i,j) += jac_ij * dFdU(1,1);
142 
143  if( this->_flow_vars.dim() == 3 )
144  {
145  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
146  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
147  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
148  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
149  (*Kww)(i,j) += jac_ij * dFdU(2,2);
150  }
151  }
152  }
153  }
154  }
155  }
156 
157 } // namespace GRINS
158 
159 // Instantiate
160 INSTANTIATE_INC_NS_SUBCLASS(VelocityDrag);
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
Definition: velocity_drag.C:46
GRINS namespace.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
Definition: velocity_drag.C:55
INSTANTIATE_INC_NS_SUBCLASS(VelocityDrag)

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