GRINS-0.6.0
parsed_velocity_source.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
28 
29 // GRINS
33 
34 // libMesh
35 #include "libmesh/quadrature.h"
36 
37 namespace GRINS
38 {
39 
40  template<class Mu>
41  ParsedVelocitySource<Mu>::ParsedVelocitySource( const std::string& physics_name, const GetPot& input )
42  : ParsedVelocitySourceBase<Mu>(physics_name, input),
43  _parsed_velocity_source_x_index(0),
44  _parsed_velocity_source_y_index(0),
45  _parsed_velocity_source_z_index(0)
46  {
47  return;
48  }
49 
50  template<class Mu>
52  {
53  return;
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  template<class Mu>
68  {
69  std::string section = "Physics/"+this->_physics_name+"/output_vars";
70 
71  std::string vel_source = "vel_source";
72 
73  if( input.have_variable(section) )
74  {
75  unsigned int n_vars = input.vector_variable_size(section);
76 
77  for( unsigned int v = 0; v < n_vars; v++ )
78  {
79  std::string name = input(section,"DIE!",v);
80 
81  if( name == std::string("parsed_velocity_source") )
82  {
83  _parsed_velocity_source_x_index =
84  postprocessing.register_quantity( vel_source+"_x" );
85 
86  _parsed_velocity_source_y_index =
87  postprocessing.register_quantity( vel_source+"_y" );
88 
89  _parsed_velocity_source_z_index =
90  postprocessing.register_quantity( vel_source+"_z" );
91  }
92  else
93  {
94  std::cerr << "Error: Invalue output_vars value for "+this->_physics_name << std::endl
95  << " Found " << name << std::endl
96  << " Acceptable values are: parsed_velocity_source" << std::endl;
97  libmesh_error();
98  }
99  }
100  }
101 
102  return;
103  }
104 
105  template<class Mu>
107  AssemblyContext& context,
108  CachedValues& /* cache */ )
109  {
110 #ifdef GRINS_USE_GRVY_TIMERS
111  this->_timer->BeginTimer("ParsedVelocitySource::element_time_derivative");
112 #endif
113 
114  // Element Jacobian * quadrature weights for interior integration
115  const std::vector<libMesh::Real> &JxW =
116  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
117 
118  // The shape functions at interior quadrature points.
119  const std::vector<std::vector<libMesh::Real> >& u_phi =
120  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
121 
122  const std::vector<libMesh::Point>& u_qpoint =
123  context.get_element_fe(this->_flow_vars.u_var())->get_xyz();
124 
125  // The number of local degrees of freedom in each variable
126  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
127 
128  // The subvectors and submatrices we need to fill:
129  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.u_var()); // R_{u},{u}
130  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.v_var()); // R_{u},{v}
131  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.u_var()); // R_{v},{u}
132  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.v_var()); // R_{v},{v}
133 
134  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
135  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
136  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
137  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
138  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
139 
140  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u_var()); // R_{u}
141  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v_var()); // R_{v}
142  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
143 
144  if( this->_dim == 3 )
145  {
146  Kuw = &context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.w_var()); // R_{u},{w}
147  Kvw = &context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.w_var()); // R_{v},{w}
148 
149  Kwu = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.u_var()); // R_{w},{u}
150  Kwv = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.v_var()); // R_{w},{v}
151  Kww = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.w_var()); // R_{w},{w}
152  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
153  }
154 
155  unsigned int n_qpoints = context.get_element_qrule().n_points();
156 
157  for (unsigned int qp=0; qp != n_qpoints; qp++)
158  {
159  libMesh::NumberVectorValue F;
160  libMesh::NumberTensorValue dFdU;
161  libMesh::NumberTensorValue* dFdU_ptr =
162  compute_jacobian ? &dFdU : NULL;
163  if (!this->compute_force(u_qpoint[qp], context.time, context, F, dFdU_ptr))
164  continue;
165 
166  const libMesh::Real jac = JxW[qp];
167 
168  for (unsigned int i=0; i != n_u_dofs; i++)
169  {
170  const libMesh::Number jac_i = jac * u_phi[i][qp];
171 
172  Fu(i) += F(0)*jac_i;
173 
174  Fv(i) += F(1)*jac_i;
175  if( this->_dim == 3 )
176  {
177  (*Fw)(i) += F(2)*jac_i;
178  }
179 
180  if( compute_jacobian )
181  {
182  for (unsigned int j=0; j != n_u_dofs; j++)
183  {
184  const libMesh::Number jac_ij = context.get_elem_solution_derivative() * jac_i * u_phi[j][qp];
185  Kuu(i,j) += jac_ij * dFdU(0,0);
186  Kuv(i,j) += jac_ij * dFdU(0,1);
187  Kvu(i,j) += jac_ij * dFdU(1,0);
188  Kvv(i,j) += jac_ij * dFdU(1,1);
189 
190  if( this->_dim == 3 )
191  {
192  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
193  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
194 
195  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
196  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
197  (*Kww)(i,j) += jac_ij * dFdU(2,2);
198  }
199  }
200  }
201  }
202  }
203 
204 
205 #ifdef GRINS_USE_GRVY_TIMERS
206  this->_timer->EndTimer("ParsedVelocitySource::element_time_derivative");
207 #endif
208 
209  return;
210  }
211 
212  template<class Mu>
214  const AssemblyContext& context,
215  const libMesh::Point& point,
216  libMesh::Real& value )
217  {
218  libMesh::DenseVector<libMesh::Number> output_vec(3);
219 
220  if( quantity_index == this->_parsed_velocity_source_x_index )
221  {
222  (*this->velocity_source_function)(context, point, context.time, output_vec);
223 
224  value = output_vec(0);
225  }
226  else if( quantity_index == this->_parsed_velocity_source_y_index )
227  {
228  (*this->velocity_source_function)(context, point, context.time, output_vec);
229 
230  value = output_vec(1);
231  }
232  else if( quantity_index == this->_parsed_velocity_source_z_index )
233  {
234  (*this->velocity_source_function)(context, point, context.time, output_vec);
235 
236  value = output_vec(2);
237  }
238 
239  return;
240  }
241 
242 } // namespace GRINS
243 // Instantiate
244 INSTANTIATE_INC_NS_SUBCLASS(ParsedVelocitySource);
virtual void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Register postprocessing variables for ParsedVelocitySource.
unsigned int register_quantity(std::string name)
Register quantity to be postprocessed.
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
GRINS namespace.
INSTANTIATE_INC_NS_SUBCLASS(ParsedVelocitySource)

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