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

Generated on Thu Jun 2 2016 21:52:27 for GRINS-0.7.0 by  doxygen 1.8.10