41 #include "libmesh/quadrature.h"
42 #include "libmesh/elem.h"
53 _spalart_allmaras_helper(input),
55 _no_of_walls(input(
"Physics/"+
PhysicsNaming::spalart_allmaras()+
"/no_of_walls", 0)),
56 _infinite_distance(input(
"Physics/"+
PhysicsNaming::spalart_allmaras()+
"/infinite_distance", false))
78 this->boundary_mesh.reset(
new libMesh::SerialMesh(system->get_mesh().comm() , system->get_mesh().mesh_dimension()) );
81 (system->get_mesh()).boundary_info->sync(_wall_ids, *boundary_mesh);
84 this->distance_function.reset(
new DistanceFunction(system->get_equation_systems(), *boundary_mesh));
98 context.get_element_fe(_turbulence_vars.nu())->get_JxW();
99 context.get_element_fe(_turbulence_vars.nu())->get_phi();
100 context.get_element_fe(_turbulence_vars.nu())->get_dphi();
101 context.get_element_fe(_turbulence_vars.nu())->get_xyz();
103 context.get_element_fe(_turbulence_vars.nu())->get_phi();
104 context.get_element_fe(_turbulence_vars.nu())->get_xyz();
106 context.get_side_fe(_turbulence_vars.nu())->get_JxW();
107 context.get_side_fe(_turbulence_vars.nu())->get_phi();
108 context.get_side_fe(_turbulence_vars.nu())->get_dphi();
109 context.get_side_fe(_turbulence_vars.nu())->get_xyz();
119 system->time_evolving(this->_turbulence_vars.nu(), 1);
128 libMesh::Elem &elem_pointer = context.get_elem();
134 const std::vector<libMesh::Real> &JxW =
135 context.get_element_fe(this->_turbulence_vars.nu())->get_JxW();
138 const std::vector<std::vector<libMesh::Real> >& nu_phi =
139 context.get_element_fe(this->_turbulence_vars.nu())->get_phi();
143 const std::vector<std::vector<libMesh::RealGradient> >& nu_gradphi =
144 context.get_element_fe(this->_turbulence_vars.nu())->get_dphi();
147 const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu()).size();
157 libMesh::DenseSubVector<libMesh::Number> &Fnu = context.get_elem_residual(this->_turbulence_vars.nu());
165 unsigned int n_qpoints = context.get_element_qrule().n_points();
168 libMesh::UniquePtr< libMesh::DenseVector<libMesh::Real> > distance_qp;
171 distance_qp = this->distance_function->interpolate(&elem_pointer, context.get_element_qrule().get_points());
173 for (
unsigned int qp=0; qp != n_qpoints; qp++)
177 nu = context.interior_value(this->_turbulence_vars.nu(), qp);
179 libMesh::Gradient grad_nu;
180 grad_nu = context.interior_gradient(this->_turbulence_vars.nu(), qp);
182 libMesh::Real jac = JxW[qp];
185 libMesh::Real mu_qp = this->_mu(context, qp);
188 libMesh::Real vorticity_value_qp = this->_spalart_allmaras_helper.vorticity(context, qp);
192 u = context.interior_value(this->_flow_vars.u(), qp);
193 v = context.interior_value(this->_flow_vars.v(), qp);
195 libMesh::NumberVectorValue U(u,v);
196 if (this->_flow_vars.dim() == 3)
197 U(2) = context.interior_value(this->_flow_vars.w(), qp);
200 libMesh::Real S_tilde = this->_sa_params.source_fn(nu, mu_qp, (*distance_qp)(qp), vorticity_value_qp, _infinite_distance);
203 libMesh::Real chi = nu/mu_qp;
204 libMesh::Real f_t2 = this->_sa_params.get_c_t3()*exp(-this->_sa_params.get_c_t4()*chi*chi);
206 libMesh::Real source_term = this->_sa_params.get_cb1()*(1 - f_t2)*S_tilde*nu;
210 source_term = this->_sa_params.get_cb1()*(1 - this->_sa_params.get_c_t3())*vorticity_value_qp*nu;
214 libMesh::Real fw = this->_sa_params.destruction_fn(nu, (*distance_qp)(qp), S_tilde, _infinite_distance);
216 libMesh::Real nud = 0.0;
217 if(_infinite_distance)
223 nud = nu/(*distance_qp)(qp);
225 libMesh::Real nud2 = nud*nud;
226 libMesh::Real kappa2 = (this->_sa_params.get_kappa())*(this->_sa_params.get_kappa());
227 libMesh::Real cw1 = this->_sa_params.get_cb1()/kappa2 + (1.0 + this->_sa_params.get_cb2())/this->_sa_params.get_sigma();
228 libMesh::Real destruction_term = (cw1*fw - (this->_sa_params.get_cb1()/kappa2)*f_t2)*nud2;
233 destruction_term = -cw1*nud2;
236 libMesh::Real fn1 = 1.0;
240 libMesh::Real chi3 = chi*chi*chi;
241 fn1 = (this->_sa_params.get_c_n1() + chi3)/(this->_sa_params.get_c_n1() - chi3);
245 for (
unsigned int i=0; i != n_nu_dofs; i++)
248 ( -this->_rho*(U*grad_nu)*nu_phi[i][qp]
249 +source_term*nu_phi[i][qp]
250 + (1./this->_sa_params.get_sigma())*(-(mu_qp+(fn1*nu))*grad_nu*nu_gradphi[i][qp] + this->_sa_params.get_cb2()*grad_nu*grad_nu*nu_phi[i][qp])
251 - destruction_term*nu_phi[i][qp]);
254 if (compute_jacobian)
256 libmesh_not_implemented();
271 const std::vector<libMesh::Real> &JxW =
272 context.get_element_fe(this->_turbulence_vars.nu())->get_JxW();
275 const std::vector<std::vector<libMesh::Real> >& nu_phi =
276 context.get_element_fe(this->_turbulence_vars.nu())->get_phi();
279 const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu()).size();
282 libMesh::DenseSubVector<libMesh::Real> &F = context.get_elem_residual(this->_turbulence_vars.nu());
286 unsigned int n_qpoints = context.get_element_qrule().n_points();
288 for (
unsigned int qp = 0; qp != n_qpoints; ++qp)
295 libMesh::Real nu_dot;
296 context.interior_rate(this->_turbulence_vars.nu(), qp, nu_dot);
298 for (
unsigned int i = 0; i != n_nu_dofs; ++i)
300 F(i) += -JxW[qp]*this->_rho*nu_dot*nu_phi[i][qp];
302 if( compute_jacobian )
304 libmesh_not_implemented();
314 (
const std::string & param_name,
319 this->_mu.register_parameter(param_name, param_pointer);
320 this->_sa_params.register_parameter(param_name, param_pointer);
virtual void init_context(AssemblyContext &context)
Initialize context for added physics variables.
GRINS::ICHandlingBase * _ic_handler
TurbulenceFEVariables & _turbulence_vars
static PhysicsName spalart_allmaras()
void check_var_subdomain_consistency(const FEVariablesBase &var) const
Check that var is enabled on at least the subdomains this Physics is.
PressureFEVariable & _press_var
Base class for reading and handling initial conditions for physics classes.
Physics class for Turbulence Models.
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.
INSTANTIATE_TURBULENCE_MODELS_SUBCLASS(SpalartAllmaras)
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
std::set< libMesh::boundary_id_type > _wall_ids
virtual void register_parameter(const std::string ¶m_name, libMesh::ParameterMultiAccessor< libMesh::Number > ¶m_pointer) const
Each subclass will register its copy of an independent.
virtual void init_variables(libMesh::FEMSystem *)
Initialize variables for this physics.
VelocityVariable & _flow_vars
void set_is_constraint_var(bool is_constraint_var)
Set whether or not this is a "constraint" variable.
virtual void mass_residual(bool compute_jacobian, AssemblyContext &context)
Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part...
unsigned int _no_of_walls
virtual void register_parameter(const std::string ¶m_name, libMesh::ParameterMultiAccessor< libMesh::Number > ¶m_pointer) const
Each subclass will register its copy of an independent.
virtual void set_time_evolving_vars(libMesh::FEMSystem *system)
Sets velocity variables to be time-evolving.