GRINS-0.8.0
velocity_penalty.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_penalty.h"
28 
29 // GRINS
33 
34 // libMesh
35 #include "libmesh/quadrature.h"
36 
37 namespace GRINS
38 {
39 
40  template<class Mu>
41  VelocityPenalty<Mu>::VelocityPenalty( const std::string& physics_name, const GetPot& input )
42  : VelocityPenaltyBase<Mu>(physics_name, input),
43  _velocity_penalty_x_index(0),
44  _velocity_penalty_y_index(0),
45  _velocity_penalty_z_index(0),
46  _velocity_penalty_base_x_index(0),
47  _velocity_penalty_base_y_index(0),
48  _velocity_penalty_base_z_index(0)
49  {
50  return;
51  }
52 
53  template<class Mu>
55  {
56  return;
57  }
58 
59  template<class Mu>
61  {
62  context.get_element_fe(this->_flow_vars.u())->get_xyz();
63  context.get_element_fe(this->_flow_vars.u())->get_phi();
64 
65  return;
66  }
67 
68  template<class Mu>
71  {
72  std::string section = "Physics/"+this->_physics_name+"/output_vars";
73 
74  std::string vel_penalty = "vel_penalty";
75  if (this->_physics_name == "VelocityPenalty2")
76  vel_penalty += '2';
77 
78  if (this->_physics_name == "VelocityPenalty3")
79  vel_penalty += '3';
80 
81  if( input.have_variable(section) )
82  {
83  unsigned int n_vars = input.vector_variable_size(section);
84 
85  for( unsigned int v = 0; v < n_vars; v++ )
86  {
87  std::string name = input(section,"DIE!",v);
88 
89  if( name == std::string("velocity_penalty") )
90  {
91  _velocity_penalty_x_index =
92  postprocessing.register_quantity( vel_penalty+"_x" );
93 
94  _velocity_penalty_y_index =
95  postprocessing.register_quantity( vel_penalty+"_y" );
96 
97  _velocity_penalty_z_index =
98  postprocessing.register_quantity( vel_penalty+"_z" );
99  }
100  else if( name == std::string("velocity_penalty_base") )
101  {
102  _velocity_penalty_base_x_index =
103  postprocessing.register_quantity( vel_penalty+"_base_x" );
104 
105  _velocity_penalty_base_y_index =
106  postprocessing.register_quantity( vel_penalty+"_base_y" );
107 
108  _velocity_penalty_base_z_index =
109  postprocessing.register_quantity( vel_penalty+"_base_z" );
110  }
111  else
112  {
113  std::cerr << "Error: Invalue output_vars value for "+this->_physics_name << std::endl
114  << " Found " << name << std::endl
115  << " Acceptable values are: velocity_penalty" << std::endl
116  << " velocity_penalty_base" << std::endl;
117  libmesh_error();
118  }
119  }
120  }
121 
122  return;
123  }
124 
125  template<class Mu>
127  ( bool compute_jacobian,
128  AssemblyContext & context )
129  {
130  // Element Jacobian * quadrature weights for interior integration
131  const std::vector<libMesh::Real> &JxW =
132  context.get_element_fe(this->_flow_vars.u())->get_JxW();
133 
134  // The shape functions at interior quadrature points.
135  const std::vector<std::vector<libMesh::Real> >& u_phi =
136  context.get_element_fe(this->_flow_vars.u())->get_phi();
137 
138  const std::vector<libMesh::Point>& u_qpoint =
139  context.get_element_fe(this->_flow_vars.u())->get_xyz();
140 
141  // The number of local degrees of freedom in each variable
142  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
143 
144  // The subvectors and submatrices we need to fill:
145  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
146  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // R_{u},{v}
147  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // R_{v},{u}
148  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
149 
150  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
151  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
152  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
153  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
154  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
155 
156  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
157  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
158  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
159 
160  if( this->_flow_vars.dim() == 3 )
161  {
162  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
163  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
164 
165  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
166  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
167  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
168  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
169  }
170 
171  unsigned int n_qpoints = context.get_element_qrule().n_points();
172 
173  for (unsigned int qp=0; qp != n_qpoints; qp++)
174  {
175  // Compute the solution at the old Newton iterate.
176  libMesh::Number u, v;
177  u = context.interior_value(this->_flow_vars.u(), qp);
178  v = context.interior_value(this->_flow_vars.v(), qp);
179 
180  libMesh::NumberVectorValue U(u,v);
181  if (this->_flow_vars.dim() == 3)
182  U(2) = context.interior_value(this->_flow_vars.w(), qp); // w
183 
184  libMesh::NumberVectorValue F;
185  libMesh::NumberTensorValue dFdU;
186  libMesh::NumberTensorValue* dFdU_ptr =
187  compute_jacobian ? &dFdU : NULL;
188  if (!this->compute_force(u_qpoint[qp], context, U, F, dFdU_ptr))
189  continue;
190 
191  const libMesh::Real jac = JxW[qp];
192 
193  for (unsigned int i=0; i != n_u_dofs; i++)
194  {
195  const libMesh::Number jac_i = jac * u_phi[i][qp];
196 
197  Fu(i) += F(0)*jac_i;
198 
199  Fv(i) += F(1)*jac_i;
200  if( this->_flow_vars.dim() == 3 )
201  {
202  (*Fw)(i) += F(2)*jac_i;
203  }
204 
205  if( compute_jacobian )
206  {
207  for (unsigned int j=0; j != n_u_dofs; j++)
208  {
209  const libMesh::Number jac_ij = context.get_elem_solution_derivative() * jac_i * u_phi[j][qp];
210  Kuu(i,j) += jac_ij * dFdU(0,0);
211  Kuv(i,j) += jac_ij * dFdU(0,1);
212  Kvu(i,j) += jac_ij * dFdU(1,0);
213  Kvv(i,j) += jac_ij * dFdU(1,1);
214 
215  if( this->_flow_vars.dim() == 3 )
216  {
217  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
218  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
219 
220  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
221  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
222  (*Kww)(i,j) += jac_ij * dFdU(2,2);
223  }
224  }
225  }
226  }
227  }
228  }
229 
230  template<class Mu>
231  void VelocityPenalty<Mu>::compute_postprocessed_quantity( unsigned int quantity_index,
232  const AssemblyContext& context,
233  const libMesh::Point& point,
234  libMesh::Real& value )
235  {
236  libMesh::DenseVector<libMesh::Number> output_vec(3);
237 
238  if( quantity_index == this->_velocity_penalty_x_index )
239  {
240  (*this->normal_vector_function)(context, point, context.time, output_vec);
241 
242  value = output_vec(0);
243  }
244  else if( quantity_index == this->_velocity_penalty_y_index )
245  {
246  (*this->normal_vector_function)(context, point, context.time, output_vec);
247 
248  value = output_vec(1);
249  }
250  else if( quantity_index == this->_velocity_penalty_z_index )
251  {
252  (*this->normal_vector_function)(context, point, context.time, output_vec);
253 
254  value = output_vec(2);
255  }
256  else if( quantity_index == this->_velocity_penalty_base_x_index )
257  {
258  (*this->base_velocity_function)(context, point, context.time, output_vec);
259 
260  value = output_vec(0);
261  }
262  else if( quantity_index == this->_velocity_penalty_base_y_index )
263  {
264  (*this->base_velocity_function)(context, point, context.time, output_vec);
265 
266  value = output_vec(1);
267  }
268  else if( quantity_index == this->_velocity_penalty_base_z_index )
269  {
270  (*this->base_velocity_function)(context, point, context.time, output_vec);
271 
272  value = output_vec(2);
273  }
274 
275  return;
276  }
277 
278 } // namespace GRINS
279 // Instantiate
280 INSTANTIATE_INC_NS_SUBCLASS(VelocityPenalty);
INSTANTIATE_INC_NS_SUBCLASS(VelocityPenalty)
unsigned int register_quantity(std::string name)
Register quantity to be postprocessed.
GRINS namespace.
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.
virtual void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Register postprocessing variables for VelocityPenalty.
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)

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