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

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