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"
67 int num_ids = input.vector_variable_size(
"QoI/Vorticity/enabled_subdomains" );
71 std::cerr <<
"Error: Must specify at least one subdomain id on which to compute vorticity." << std::endl;
75 for(
int i = 0; i < num_ids; i++ )
77 libMesh::subdomain_id_type s_id = input(
"QoI/Vorticity/enabled_subdomains", -1, i );
78 _subdomain_ids.insert( s_id );
86 libMesh::FEBase* u_fe = NULL;
87 libMesh::FEBase* v_fe = NULL;
89 context.get_element_fe<libMesh::Real>(
_flow_vars->
u(), u_fe);
90 context.get_element_fe<libMesh::Real>(
_flow_vars->
v(), v_fe);
101 const unsigned int qoi_index )
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();
110 unsigned int n_qpoints = context.get_element_qrule().n_points();
113 libMesh::Number& qoi = context.get_qois()[qoi_index];
115 for(
unsigned int qp = 0; qp != n_qpoints; qp++ )
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];
129 const unsigned int qoi_index )
135 libMesh::FEBase* element_fe;
136 context.get_element_fe<libMesh::Real>(
_flow_vars->
u(), element_fe);
139 const std::vector<libMesh::Real> &JxW = element_fe->get_JxW();
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();
148 const unsigned int n_T_dofs = context.get_dof_indices(0).size();
149 unsigned int n_qpoints = context.get_element_qrule().n_points();
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());
158 for(
unsigned int qp = 0; qp != n_qpoints; qp++ )
160 for(
unsigned int i = 0; i != n_T_dofs; i++ )
162 Qu(i) += - dv_phi[i][qp](1) * JxW[qp];
163 Qv(i) += du_phi[i][qp](0) * JxW[qp];
Vorticity()
User never call default constructor.
virtual QoIBase * clone() const
Required to provide clone (deep-copy) for adding QoI object to libMesh objects.
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)
std::set< libMesh::subdomain_id_type > _subdomain_ids
List of sumdomain ids for which we want to compute this QoI.
const VelocityVariable * _flow_vars
virtual void element_qoi(AssemblyContext &context, const unsigned int qoi_index)
Compute the qoi value.
Interface with libMesh for solving Multiphysics problems.
virtual void init(const GetPot &input, const MultiphysicsSystem &system, unsigned int qoi_num)
Initialize local variables.
virtual void element_qoi_derivative(AssemblyContext &context, const unsigned int qoi_index)
Compute the qoi derivative with respect to the solution.