GRINS-0.7.0
vorticity.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/vorticity.h"
28 
29 // GRINS
30 #include "grins/multiphysics_sys.h"
31 #include "grins/assembly_context.h"
32 
33 // libMesh
34 #include "libmesh/getpot.h"
35 #include "libmesh/fem_system.h"
36 #include "libmesh/quadrature.h"
37 #include "libmesh/fe_base.h"
38 #include "libmesh/elem.h"
39 
40 namespace GRINS
41 {
42  Vorticity::Vorticity( const std::string& qoi_name )
43  : QoIBase(qoi_name)
44  {
45  return;
46  }
47 
49  {
50  return;
51  }
52 
54  {
55  return new Vorticity( *this );
56  }
57 
58  void Vorticity::init
59  (const GetPot& input,
60  const MultiphysicsSystem& system,
61  unsigned int /*qoi_num*/ )
62  {
63  // Extract subdomain on which to compute to qoi
64  int num_ids = input.vector_variable_size( "QoI/Vorticity/enabled_subdomains" );
65 
66  if( num_ids == 0 )
67  {
68  std::cerr << "Error: Must specify at least one subdomain id on which to compute vorticity." << std::endl;
69  libmesh_error();
70  }
71 
72  for( int i = 0; i < num_ids; i++ )
73  {
74  libMesh::subdomain_id_type s_id = input( "QoI/Vorticity/enabled_subdomains", -1, i );
75  _subdomain_ids.insert( s_id );
76  }
77 
78  // Grab velocity variable indices
79  std::string u_var_name = input("Physics/VariableNames/u_velocity", u_var_name_default);
80  std::string v_var_name = input("Physics/VariableNames/v_velocity", v_var_name_default);
81  this->_u_var = system.variable_number(u_var_name);
82  this->_v_var = system.variable_number(v_var_name);
83 
84  return;
85  }
86 
88  {
89  libMesh::FEBase* u_fe = NULL;
90  libMesh::FEBase* v_fe = NULL;
91 
92  context.get_element_fe<libMesh::Real>(this->_u_var, u_fe);
93  context.get_element_fe<libMesh::Real>(this->_v_var, v_fe);
94 
95  u_fe->get_dphi();
96  u_fe->get_JxW();
97 
98  v_fe->get_dphi();
99 
100  return;
101  }
102 
104  const unsigned int qoi_index )
105  {
106 
107  if( _subdomain_ids.find( (&context.get_elem())->subdomain_id() ) != _subdomain_ids.end() )
108  {
109  libMesh::FEBase* element_fe;
110  context.get_element_fe<libMesh::Real>(this->_u_var, element_fe);
111  const std::vector<libMesh::Real> &JxW = element_fe->get_JxW();
112 
113  unsigned int n_qpoints = context.get_element_qrule().n_points();
114 
116  libMesh::Number& qoi = context.get_qois()[qoi_index];
117 
118  for( unsigned int qp = 0; qp != n_qpoints; qp++ )
119  {
120  libMesh::Gradient grad_u = 0.;
121  libMesh::Gradient grad_v = 0.;
122  context.interior_gradient( this->_u_var, qp, grad_u );
123  context.interior_gradient( this->_v_var, qp, grad_v );
124  qoi += (grad_v(0) - grad_u(1)) * JxW[qp];
125  }
126  }
127 
128  return;
129  }
130 
132  const unsigned int qoi_index )
133  {
134 
135  if( _subdomain_ids.find( (&context.get_elem())->subdomain_id() ) != _subdomain_ids.end() )
136  {
137  // Element
138  libMesh::FEBase* element_fe;
139  context.get_element_fe<libMesh::Real>(this->_u_var, element_fe);
140 
141  // Jacobian times integration weights
142  const std::vector<libMesh::Real> &JxW = element_fe->get_JxW();
143 
144  // Grad of basis functions
145  const std::vector<std::vector<libMesh::RealGradient> >& du_phi =
146  context.get_element_fe(_u_var)->get_dphi();
147  const std::vector<std::vector<libMesh::RealGradient> >& dv_phi =
148  context.get_element_fe(_v_var)->get_dphi();
149 
150  // Local DOF count and quadrature point count
151  const unsigned int n_T_dofs = context.get_dof_indices(0).size();
152  unsigned int n_qpoints = context.get_element_qrule().n_points();
153 
154  // Warning: we assume here that vorticity is the only QoI!
155  // This should be consistent with the assertion in grins_mesh_adaptive_solver.C
157  libMesh::DenseSubVector<libMesh::Number> &Qu = context.get_qoi_derivatives(qoi_index, _u_var);
158  libMesh::DenseSubVector<libMesh::Number> &Qv = context.get_qoi_derivatives(qoi_index, _v_var);
159 
160  // Integration loop
161  for( unsigned int qp = 0; qp != n_qpoints; qp++ )
162  {
163  for( unsigned int i = 0; i != n_T_dofs; i++ )
164  {
165  Qu(i) += - dv_phi[i][qp](1) * JxW[qp];
166  Qv(i) += du_phi[i][qp](0) * JxW[qp];
167  }
168  }
169  }
170 
171  return;
172  }
173 
174 } //namespace GRINS
Vorticity()
User never call default constructor.
virtual ~Vorticity()
Definition: vorticity.C:48
virtual QoIBase * clone() const
Required to provide clone (deep-copy) for adding QoI object to libMesh objects.
Definition: vorticity.C:53
virtual void init_context(AssemblyContext &context)
Definition: vorticity.C:87
std::set< libMesh::subdomain_id_type > _subdomain_ids
List of sumdomain ids for which we want to compute this QoI.
Definition: vorticity.h:91
const std::string v_var_name_default
y-velocity
GRINS namespace.
const std::string u_var_name_default
Default physics variable names.
virtual void element_qoi(AssemblyContext &context, const unsigned int qoi_index)
Compute the qoi value.
Definition: vorticity.C:103
Interface with libMesh for solving Multiphysics problems.
virtual void init(const GetPot &input, const MultiphysicsSystem &system, unsigned int qoi_num)
Initialize local variables.
Definition: vorticity.C:59
virtual void element_qoi_derivative(AssemblyContext &context, const unsigned int qoi_index)
Compute the qoi derivative with respect to the solution.
Definition: vorticity.C:131
VariableIndex _v_var
v-velocity component variable index
Definition: vorticity.h:88
VariableIndex _u_var
u-velocity component variable index
Definition: vorticity.h:85

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