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

Generated on Tue Dec 19 2017 12:47:28 for GRINS-0.8.0 by  doxygen 1.8.9.1