GRINS-0.6.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-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/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  this->read_input_options(input);
45 
46  return;
47  }
48 
49  template<class Mu>
51  {
52  return;
53  }
54 
55 
56  template<class Mu>
58  {
59  context.get_element_fe(this->_flow_vars.u_var())->get_xyz();
60  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
61 
62  return;
63  }
64 
65 
66  template<class Mu>
67  void VelocityDrag<Mu>::element_time_derivative( bool compute_jacobian,
68  AssemblyContext& context,
69  CachedValues& /* cache */ )
70  {
71 #ifdef GRINS_USE_GRVY_TIMERS
72  this->_timer->BeginTimer("VelocityDrag::element_time_derivative");
73 #endif
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_var())->get_JxW();
78 
79  // The shape functions at interior quadrature points.
80  const std::vector<std::vector<libMesh::Real> >& u_phi =
81  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
82 
83  const std::vector<libMesh::Point>& u_qpoint =
84  context.get_element_fe(this->_flow_vars.u_var())->get_xyz();
85 
86  // The number of local degrees of freedom in each variable
87  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
88 
89  // The subvectors and submatrices we need to fill:
90  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.u_var()); // R_{u},{u}
91  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.v_var()); // R_{u},{v}
92  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.u_var()); // R_{v},{u}
93  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.v_var()); // R_{v},{v}
94 
95  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
96  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
97  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
98  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
99  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
100 
101  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u_var()); // R_{u}
102  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v_var()); // R_{v}
103  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
104 
105  if( this->_dim == 3 )
106  {
107  Kuw = &context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.w_var()); // R_{u},{w}
108  Kvw = &context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.w_var()); // R_{v},{w}
109 
110  Kwu = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.u_var()); // R_{w},{u}
111  Kwv = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.v_var()); // R_{w},{v}
112  Kww = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.w_var()); // R_{w},{w}
113  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
114  }
115 
116  unsigned int n_qpoints = context.get_element_qrule().n_points();
117 
118  for (unsigned int qp=0; qp != n_qpoints; qp++)
119  {
120  // Compute the solution at the old Newton iterate.
121  libMesh::Number u, v;
122  u = context.interior_value(this->_flow_vars.u_var(), qp);
123  v = context.interior_value(this->_flow_vars.v_var(), qp);
124 
125  libMesh::NumberVectorValue U(u,v);
126  if (this->_dim == 3)
127  U(2) = context.interior_value(this->_flow_vars.w_var(), qp); // w
128 
129  libMesh::NumberVectorValue F;
130  libMesh::NumberTensorValue dFdU;
131  libMesh::NumberTensorValue* dFdU_ptr =
132  compute_jacobian ? &dFdU : NULL;
133  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
134  continue;
135 
136  libMesh::Real jac = JxW[qp];
137 
138  for (unsigned int i=0; i != n_u_dofs; i++)
139  {
140  const libMesh::Number jac_i = jac * u_phi[i][qp];
141 
142  Fu(i) += F(0)*jac_i;
143  Fv(i) += F(1)*jac_i;
144 
145  if( this->_dim == 3 )
146  (*Fw)(i) += F(2)*jac_i;
147 
148  if( compute_jacobian )
149  {
150  for (unsigned int j=0; j != n_u_dofs; j++)
151  {
152  const libMesh::Number jac_ij =
153  jac_i * context.get_elem_solution_derivative() *
154  u_phi[j][qp];
155  Kuu(i,j) += jac_ij * dFdU(0,0);
156  Kuv(i,j) += jac_ij * dFdU(0,1);
157  Kvu(i,j) += jac_ij * dFdU(1,0);
158  Kvv(i,j) += jac_ij * dFdU(1,1);
159 
160  if( this->_dim == 3 )
161  {
162  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
163  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
164  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
165  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
166  (*Kww)(i,j) += jac_ij * dFdU(2,2);
167  }
168  }
169  }
170  }
171  }
172 
173 
174 #ifdef GRINS_USE_GRVY_TIMERS
175  this->_timer->EndTimer("VelocityDrag::element_time_derivative");
176 #endif
177 
178  return;
179  }
180 
181 } // namespace GRINS
182 
183 // Instantiate
184 INSTANTIATE_INC_NS_SUBCLASS(VelocityDrag);
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
Definition: velocity_drag.C:57
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:67
GRINS namespace.
virtual void read_input_options(const GetPot &input)
Read options from GetPot input file.
INSTANTIATE_INC_NS_SUBCLASS(VelocityDrag)

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