GRINS-0.8.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-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
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())->get_xyz();
60  context.get_element_fe(this->_flow_vars.u())->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  ( bool compute_jacobian,
108  AssemblyContext & context )
109  {
110  // Element Jacobian * quadrature weights for interior integration
111  const std::vector<libMesh::Real> &JxW =
112  context.get_element_fe(this->_flow_vars.u())->get_JxW();
113 
114  // The shape functions at interior quadrature points.
115  const std::vector<std::vector<libMesh::Real> >& u_phi =
116  context.get_element_fe(this->_flow_vars.u())->get_phi();
117 
118  const std::vector<libMesh::Point>& u_qpoint =
119  context.get_element_fe(this->_flow_vars.u())->get_xyz();
120 
121  // The number of local degrees of freedom in each variable
122  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u()).size();
123 
124  // The subvectors and submatrices we need to fill:
125  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.u()); // R_{u},{u}
126  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.v()); // R_{u},{v}
127  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.u()); // R_{v},{u}
128  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.v()); // R_{v},{v}
129 
130  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
131  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
132  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
133  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
134  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
135 
136  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u()); // R_{u}
137  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v()); // R_{v}
138  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
139 
140  if( this->_flow_vars.dim() == 3 )
141  {
142  Kuw = &context.get_elem_jacobian(this->_flow_vars.u(), this->_flow_vars.w()); // R_{u},{w}
143  Kvw = &context.get_elem_jacobian(this->_flow_vars.v(), this->_flow_vars.w()); // R_{v},{w}
144 
145  Kwu = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.u()); // R_{w},{u}
146  Kwv = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.v()); // R_{w},{v}
147  Kww = &context.get_elem_jacobian(this->_flow_vars.w(), this->_flow_vars.w()); // R_{w},{w}
148  Fw = &context.get_elem_residual(this->_flow_vars.w()); // R_{w}
149  }
150 
151  unsigned int n_qpoints = context.get_element_qrule().n_points();
152 
153  for (unsigned int qp=0; qp != n_qpoints; qp++)
154  {
155  libMesh::NumberVectorValue F;
156  libMesh::NumberTensorValue dFdU;
157  libMesh::NumberTensorValue* dFdU_ptr =
158  compute_jacobian ? &dFdU : NULL;
159  if (!this->compute_force(u_qpoint[qp], context.time, context, F, dFdU_ptr))
160  continue;
161 
162  const libMesh::Real jac = JxW[qp];
163 
164  for (unsigned int i=0; i != n_u_dofs; i++)
165  {
166  const libMesh::Number jac_i = jac * u_phi[i][qp];
167 
168  Fu(i) += F(0)*jac_i;
169 
170  Fv(i) += F(1)*jac_i;
171  if( this->_flow_vars.dim() == 3 )
172  {
173  (*Fw)(i) += F(2)*jac_i;
174  }
175 
176  if( compute_jacobian )
177  {
178  for (unsigned int j=0; j != n_u_dofs; j++)
179  {
180  const libMesh::Number jac_ij = context.get_elem_solution_derivative() * jac_i * u_phi[j][qp];
181  Kuu(i,j) += jac_ij * dFdU(0,0);
182  Kuv(i,j) += jac_ij * dFdU(0,1);
183  Kvu(i,j) += jac_ij * dFdU(1,0);
184  Kvv(i,j) += jac_ij * dFdU(1,1);
185 
186  if( this->_flow_vars.dim() == 3 )
187  {
188  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
189  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
190 
191  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
192  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
193  (*Kww)(i,j) += jac_ij * dFdU(2,2);
194  }
195  }
196  }
197  }
198  }
199  }
200 
201  template<class Mu>
203  const AssemblyContext& context,
204  const libMesh::Point& point,
205  libMesh::Real& value )
206  {
207  libMesh::DenseVector<libMesh::Number> output_vec(3);
208 
209  if( quantity_index == this->_parsed_velocity_source_x_index )
210  {
211  (*this->velocity_source_function)(context, point, context.time, output_vec);
212 
213  value = output_vec(0);
214  }
215  else if( quantity_index == this->_parsed_velocity_source_y_index )
216  {
217  (*this->velocity_source_function)(context, point, context.time, output_vec);
218 
219  value = output_vec(1);
220  }
221  else if( quantity_index == this->_parsed_velocity_source_z_index )
222  {
223  (*this->velocity_source_function)(context, point, context.time, output_vec);
224 
225  value = output_vec(2);
226  }
227 
228  return;
229  }
230 
231 } // namespace GRINS
232 // Instantiate
233 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)
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 Tue Dec 19 2017 12:47:28 for GRINS-0.8.0 by  doxygen 1.8.9.1