GRINS-0.6.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-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_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_var())->get_xyz();
63  context.get_element_fe(this->_flow_vars.u_var())->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>
126  void VelocityPenalty<Mu>::element_time_derivative( bool compute_jacobian,
127  AssemblyContext& context,
128  CachedValues& /* cache */ )
129  {
130 #ifdef GRINS_USE_GRVY_TIMERS
131  this->_timer->BeginTimer("VelocityPenalty::element_time_derivative");
132 #endif
133 
134  // Element Jacobian * quadrature weights for interior integration
135  const std::vector<libMesh::Real> &JxW =
136  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
137 
138  // The shape functions at interior quadrature points.
139  const std::vector<std::vector<libMesh::Real> >& u_phi =
140  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
141 
142  const std::vector<libMesh::Point>& u_qpoint =
143  context.get_element_fe(this->_flow_vars.u_var())->get_xyz();
144 
145  // The number of local degrees of freedom in each variable
146  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
147 
148  // The subvectors and submatrices we need to fill:
149  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.u_var()); // R_{u},{u}
150  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.v_var()); // R_{u},{v}
151  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.u_var()); // R_{v},{u}
152  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.v_var()); // R_{v},{v}
153 
154  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
155  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
156  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
157  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
158  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
159 
160  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u_var()); // R_{u}
161  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v_var()); // R_{v}
162  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
163 
164  if( this->_dim == 3 )
165  {
166  Kuw = &context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.w_var()); // R_{u},{w}
167  Kvw = &context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.w_var()); // R_{v},{w}
168 
169  Kwu = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.u_var()); // R_{w},{u}
170  Kwv = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.v_var()); // R_{w},{v}
171  Kww = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.w_var()); // R_{w},{w}
172  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
173  }
174 
175  unsigned int n_qpoints = context.get_element_qrule().n_points();
176 
177  for (unsigned int qp=0; qp != n_qpoints; qp++)
178  {
179  // Compute the solution at the old Newton iterate.
180  libMesh::Number u, v;
181  u = context.interior_value(this->_flow_vars.u_var(), qp);
182  v = context.interior_value(this->_flow_vars.v_var(), qp);
183 
184  libMesh::NumberVectorValue U(u,v);
185  if (this->_dim == 3)
186  U(2) = context.interior_value(this->_flow_vars.w_var(), qp); // w
187 
188  libMesh::NumberVectorValue F;
189  libMesh::NumberTensorValue dFdU;
190  libMesh::NumberTensorValue* dFdU_ptr =
191  compute_jacobian ? &dFdU : NULL;
192  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
193  continue;
194 
195  const libMesh::Real jac = JxW[qp];
196 
197  for (unsigned int i=0; i != n_u_dofs; i++)
198  {
199  const libMesh::Number jac_i = jac * u_phi[i][qp];
200 
201  Fu(i) += F(0)*jac_i;
202 
203  Fv(i) += F(1)*jac_i;
204  if( this->_dim == 3 )
205  {
206  (*Fw)(i) += F(2)*jac_i;
207  }
208 
209  if( compute_jacobian )
210  {
211  for (unsigned int j=0; j != n_u_dofs; j++)
212  {
213  const libMesh::Number jac_ij = context.get_elem_solution_derivative() * jac_i * u_phi[j][qp];
214  Kuu(i,j) += jac_ij * dFdU(0,0);
215  Kuv(i,j) += jac_ij * dFdU(0,1);
216  Kvu(i,j) += jac_ij * dFdU(1,0);
217  Kvv(i,j) += jac_ij * dFdU(1,1);
218 
219  if( this->_dim == 3 )
220  {
221  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
222  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
223 
224  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
225  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
226  (*Kww)(i,j) += jac_ij * dFdU(2,2);
227  }
228  }
229  }
230  }
231  }
232 
233 
234 #ifdef GRINS_USE_GRVY_TIMERS
235  this->_timer->EndTimer("VelocityPenalty::element_time_derivative");
236 #endif
237 
238  return;
239  }
240 
241  template<class Mu>
242  void VelocityPenalty<Mu>::compute_postprocessed_quantity( unsigned int quantity_index,
243  const AssemblyContext& context,
244  const libMesh::Point& point,
245  libMesh::Real& value )
246  {
247  libMesh::DenseVector<libMesh::Number> output_vec(3);
248 
249  if( quantity_index == this->_velocity_penalty_x_index )
250  {
251  (*this->normal_vector_function)(point, context.time, output_vec);
252 
253  value = output_vec(0);
254  }
255  else if( quantity_index == this->_velocity_penalty_y_index )
256  {
257  (*this->normal_vector_function)(point, context.time, output_vec);
258 
259  value = output_vec(1);
260  }
261  else if( quantity_index == this->_velocity_penalty_z_index )
262  {
263  (*this->normal_vector_function)(point, context.time, output_vec);
264 
265  value = output_vec(2);
266  }
267  else if( quantity_index == this->_velocity_penalty_base_x_index )
268  {
269  (*this->base_velocity_function)(point, context.time, output_vec);
270 
271  value = output_vec(0);
272  }
273  else if( quantity_index == this->_velocity_penalty_base_y_index )
274  {
275  (*this->base_velocity_function)(point, context.time, output_vec);
276 
277  value = output_vec(1);
278  }
279  else if( quantity_index == this->_velocity_penalty_base_z_index )
280  {
281  (*this->base_velocity_function)(point, context.time, output_vec);
282 
283  value = output_vec(2);
284  }
285 
286  return;
287  }
288 
289 } // namespace GRINS
290 // Instantiate
291 INSTANTIATE_INC_NS_SUBCLASS(VelocityPenalty);
INSTANTIATE_INC_NS_SUBCLASS(VelocityPenalty)
unsigned int register_quantity(std::string name)
Register quantity to be postprocessed.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
GRINS namespace.
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 Mon Jun 22 2015 21:32:20 for GRINS-0.6.0 by  doxygen 1.8.9.1