GRINS-0.7.0
List of all members | Public Member Functions | Public Attributes | Protected Attributes | Private Member Functions
GRINS::SpalartAllmaras< Viscosity > Class Template Reference

Physics class for Incompressible Navier-Stokes. More...

#include <spalart_allmaras.h>

Inheritance diagram for GRINS::SpalartAllmaras< Viscosity >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::SpalartAllmaras< Viscosity >:
Collaboration graph
[legend]

Public Member Functions

 SpalartAllmaras (const std::string &physics_name, const GetPot &input)
 
 ~SpalartAllmaras ()
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialize variables for this physics. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
- Public Member Functions inherited from GRINS::TurbulenceModelsBase< Viscosity >
 TurbulenceModelsBase (const std::string &physics_name, const GetPot &input)
 
 ~TurbulenceModelsBase ()
 
- Public Member Functions inherited from GRINS::Physics
 Physics (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~Physics ()
 
virtual bool enabled_on_elem (const libMesh::Elem *elem)
 Find if current physics is active on supplied element. More...
 
void set_is_steady (bool is_steady)
 Sets whether this physics is to be solved with a steady solver or not. More...
 
bool is_steady () const
 Returns whether or not this physics is being solved with a steady solver. More...
 
virtual void auxiliary_init (MultiphysicsSystem &system)
 Any auxillary initialization a Physics class may need. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. More...
 
virtual void side_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_damping_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 

Public Attributes

libMesh::UniquePtr< DistanceFunctiondistance_function
 
libMesh::UniquePtr< libMesh::SerialMesh > boundary_mesh
 

Protected Attributes

VelocityFEVariables _flow_vars
 
PressureFEVariable _press_var
 
TurbulenceFEVariables _turbulence_vars
 
SpalartAllmarasHelper _spalart_allmaras_helper
 
SpalartAllmarasParameters _sa_params
 Object handling the plethora of parameters. More...
 
std::set< libMesh::boundary_id_type > _wall_ids
 
unsigned int _no_of_walls
 
bool _infinite_distance
 
- Protected Attributes inherited from GRINS::TurbulenceModelsBase< Viscosity >
unsigned int _dim
 Physical dimension of problem. More...
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
Viscosity _mu
 Viscosity object. More...
 
- Protected Attributes inherited from GRINS::Physics
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 

Private Member Functions

 SpalartAllmaras ()
 
void register_variables ()
 

Additional Inherited Members

- Static Public Member Functions inherited from GRINS::Physics
static void set_is_axisymmetric (bool is_axisymmetric)
 Set whether we should treat the problem as axisymmetric. More...
 
static bool is_axisymmetric ()
 
- Static Public Attributes inherited from GRINS::ParameterUser
static std::string zero_vector_function = std::string("{0}")
 A parseable function string with LIBMESH_DIM components, all 0. More...
 
- Protected Member Functions inherited from GRINS::Physics
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe (const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
 
void parse_enabled_subdomains (const GetPot &input, const std::string &physics_name)
 
- Static Protected Attributes inherited from GRINS::Physics
static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 
static bool _is_axisymmetric = false
 Caches whether we are solving an axisymmetric problem or not. More...
 

Detailed Description

template<class Viscosity>
class GRINS::SpalartAllmaras< Viscosity >

Physics class for Incompressible Navier-Stokes.

This physics class implements the classical Incompressible Navier-Stokes equations. This is a templated class, the class Viscosity can be instantiated as a specific type (right now:ConstantViscosity or SpatiallyVaryingViscosity) to allow the user to specify a constant or spatially varying viscosity in the input file

Definition at line 58 of file spalart_allmaras.h.

Constructor & Destructor Documentation

template<class Mu >
GRINS::SpalartAllmaras< Mu >::SpalartAllmaras ( const std::string &  physics_name,
const GetPot &  input 
)

Definition at line 48 of file spalart_allmaras.C.

References GRINS::Physics::_ic_handler, GRINS::SpalartAllmaras< Viscosity >::_no_of_walls, GRINS::SpalartAllmaras< Viscosity >::_wall_ids, GRINS::SpalartAllmaras< Viscosity >::register_variables(), and GRINS::PhysicsNaming::spalart_allmaras().

49  : TurbulenceModelsBase<Mu>(physics_name, input), // Define class variables
51  _press_var(input,PhysicsNaming::incompressible_navier_stokes(), true /*is_constraint_var*/),
54  _sa_params(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))
57  {
58  // Loop over the _no_of_walls and fill the wall_ids set
59  for(unsigned int i = 0; i != _no_of_walls; i++)
60  {
61  _wall_ids.insert(input("Physics/"+PhysicsNaming::spalart_allmaras()+"/wall_ids", 0, i ));
62  }
63 
64  std::cout<<"No of walls: "<<_no_of_walls<<std::endl;
65 
66  for( std::set<libMesh::boundary_id_type>::iterator b_id = _wall_ids.begin(); b_id != _wall_ids.end(); ++b_id )
67  {
68  std::cout<<"Boundary Id: "<<*b_id<<std::endl;
69  }
70 
71  this->register_variables();
72 
73  this->_ic_handler = new GenericICHandler( physics_name, input );
74  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:269
static PhysicsName spalart_allmaras()
PressureFEVariable _press_var
SpalartAllmarasParameters _sa_params
Object handling the plethora of parameters.
SpalartAllmarasHelper _spalart_allmaras_helper
TurbulenceFEVariables _turbulence_vars
std::set< libMesh::boundary_id_type > _wall_ids
static PhysicsName incompressible_navier_stokes()
VelocityFEVariables _flow_vars
template<class Mu >
GRINS::SpalartAllmaras< Mu >::~SpalartAllmaras ( )

Definition at line 77 of file spalart_allmaras.C.

78  {
79  return;
80  }
template<class Viscosity >
GRINS::SpalartAllmaras< Viscosity >::SpalartAllmaras ( )
private

Member Function Documentation

template<class Mu >
void GRINS::SpalartAllmaras< Mu >::element_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Time dependent part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Reimplemented in GRINS::SpalartAllmarasSPGSMStabilization< Viscosity >.

Definition at line 156 of file spalart_allmaras.C.

159  {
160 #ifdef GRINS_USE_GRVY_TIMERS
161  this->_timer->BeginTimer("SpalartAllmaras::element_time_derivative");
162 #endif
163 
164  // Get a pointer to the current element, we need this for computing
165  // the distance to wall for the quadrature points
166  libMesh::Elem &elem_pointer = context.get_elem();
167 
168  // We get some references to cell-specific data that
169  // will be used to assemble the linear system.
170 
171  // Element Jacobian * quadrature weights for interior integration.
172  const std::vector<libMesh::Real> &JxW =
173  context.get_element_fe(this->_turbulence_vars.nu())->get_JxW();
174 
175  // The viscosity shape functions at interior quadrature points.
176  const std::vector<std::vector<libMesh::Real> >& nu_phi =
177  context.get_element_fe(this->_turbulence_vars.nu())->get_phi();
178 
179  // The viscosity shape function gradients (in global coords.)
180  // at interior quadrature points.
181  const std::vector<std::vector<libMesh::RealGradient> >& nu_gradphi =
182  context.get_element_fe(this->_turbulence_vars.nu())->get_dphi();
183 
184  // The number of local degrees of freedom in each variable.
185  const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu()).size();
186 
187  // The subvectors and submatrices we need to fill:
188  //
189  // K_{\alpha \beta} = R_{\alpha},{\beta} = \partial{ R_{\alpha} } / \partial{ {\beta} } (where R denotes residual)
190  // e.g., for \alpha = v and \beta = u we get: K{vu} = R_{v},{u}
191  // Note that Kpu, Kpv, Kpw and Fp comes as constraint.
192 
193  //libMesh::DenseSubMatrix<libMesh::Number> &Knunu = context.get_elem_jacobian(this->_turbulence_vars.nu(), this->_turbulence_vars.nu()); // R_{nu},{nu}
194 
195  libMesh::DenseSubVector<libMesh::Number> &Fnu = context.get_elem_residual(this->_turbulence_vars.nu()); // R_{nu}
196 
197  // Now we will build the element Jacobian and residual.
198  // Constructing the residual requires the solution and its
199  // gradient from the previous timestep. This must be
200  // calculated at each quadrature point by summing the
201  // solution degree-of-freedom values by the appropriate
202  // weight functions.
203  unsigned int n_qpoints = context.get_element_qrule().n_points();
204 
205  // Auto pointer to distance fcn evaluated at quad points
206  libMesh::UniquePtr< libMesh::DenseVector<libMesh::Real> > distance_qp;
207 
208  // Fill the vector of distances to quadrature points
209  distance_qp = this->distance_function->interpolate(&elem_pointer, context.get_element_qrule().get_points());
210 
211  for (unsigned int qp=0; qp != n_qpoints; qp++)
212  {
213  // Compute the solution & its gradient at the old Newton iterate.
214  libMesh::Number nu;
215  nu = context.interior_value(this->_turbulence_vars.nu(), qp);
216 
217  libMesh::Gradient grad_nu;
218  grad_nu = context.interior_gradient(this->_turbulence_vars.nu(), qp);
219 
220  libMesh::Real jac = JxW[qp];
221 
222  // The physical viscosity
223  libMesh::Real mu_qp = this->_mu(context, qp);
224 
225  // The vorticity value
226  libMesh::Real vorticity_value_qp = this->_spalart_allmaras_helper.vorticity(context, qp);
227 
228  // The flow velocity
229  libMesh::Number u,v;
230  u = context.interior_value(this->_flow_vars.u(), qp);
231  v = context.interior_value(this->_flow_vars.v(), qp);
232 
233  libMesh::NumberVectorValue U(u,v);
234  if (this->_dim == 3)
235  U(2) = context.interior_value(this->_flow_vars.w(), qp);
236 
237  //The source term
238  libMesh::Real S_tilde = this->_sa_params.source_fn(nu, mu_qp, (*distance_qp)(qp), vorticity_value_qp, _infinite_distance);
239 
240  // The ft2 function needed for the negative S-A model
241  libMesh::Real chi = nu/mu_qp;
242  libMesh::Real f_t2 = this->_sa_params.get_c_t3()*exp(-this->_sa_params.get_c_t4()*chi*chi);
243 
244  libMesh::Real source_term = this->_sa_params.get_cb1()*(1 - f_t2)*S_tilde*nu;
245  // For a negative turbulent viscosity nu < 0.0 we need to use a different production function
246  if(nu < 0.0)
247  {
248  source_term = this->_sa_params.get_cb1()*(1 - this->_sa_params.get_c_t3())*vorticity_value_qp*nu;
249  }
250 
251  // The wall destruction term
252  libMesh::Real fw = this->_sa_params.destruction_fn(nu, (*distance_qp)(qp), S_tilde, _infinite_distance);
253 
254  libMesh::Real nud = 0.0;
256  {
257  nud = 0.0;
258  }
259  else
260  {
261  nud = nu/(*distance_qp)(qp);
262  }
263  libMesh::Real nud2 = nud*nud;
264  libMesh::Real kappa2 = (this->_sa_params.get_kappa())*(this->_sa_params.get_kappa());
265  libMesh::Real cw1 = this->_sa_params.get_cb1()/kappa2 + (1.0 + this->_sa_params.get_cb2())/this->_sa_params.get_sigma();
266  libMesh::Real destruction_term = (cw1*fw - (this->_sa_params.get_cb1()/kappa2)*f_t2)*nud2;
267 
268  // For a negative turbulent viscosity nu < 0.0 we need to use a different production function
269  if(nu < 0.0)
270  {
271  destruction_term = -cw1*nud2;
272  }
273 
274  libMesh::Real fn1 = 1.0;
275  // For a negative turbulent viscosity, fn1 needs to be calculated
276  if(nu < 0.0)
277  {
278  libMesh::Real chi3 = chi*chi*chi;
279  fn1 = (this->_sa_params.get_c_n1() + chi3)/(this->_sa_params.get_c_n1() - chi3);
280  }
281 
282  // First, an i-loop over the viscosity degrees of freedom.
283  for (unsigned int i=0; i != n_nu_dofs; i++)
284  {
285  Fnu(i) += jac *
286  ( -this->_rho*(U*grad_nu)*nu_phi[i][qp] // convection term (assumes incompressibility)
287  +source_term*nu_phi[i][qp] // source term
288  + (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]) // diffusion term
289  - destruction_term*nu_phi[i][qp]); // destruction term
290 
291  // Compute the jacobian if not using numerical jacobians
292  if (compute_jacobian)
293  {
294  libmesh_not_implemented();
295  } // end - if (compute_jacobian)
296 
297  } // end of the outer dof (i) loop
298  } // end of the quadrature point (qp) loop
299 
300 #ifdef GRINS_USE_GRVY_TIMERS
301  this->_timer->EndTimer("SpalartAllmaras::element_time_derivative");
302 #endif
303 
304  return;
305  }
libMesh::UniquePtr< DistanceFunction > distance_function
unsigned int _dim
Physical dimension of problem.
libMesh::Real source_fn(libMesh::Number nu, libMesh::Real mu, libMesh::Real wall_distance, libMesh::Real vorticity_value, bool infinite_distance) const
SpalartAllmarasParameters _sa_params
Object handling the plethora of parameters.
Viscosity _mu
Viscosity object.
SpalartAllmarasHelper _spalart_allmaras_helper
TurbulenceFEVariables _turbulence_vars
libMesh::Real vorticity(AssemblyContext &context, unsigned int qp) const
VelocityFEVariables _flow_vars
libMesh::Real destruction_fn(libMesh::Number nu, libMesh::Real wall_distance, libMesh::Real S_tilde, bool infinite_distance) const
libMesh::Number _rho
Material parameters, read from input.
template<class Mu >
void GRINS::SpalartAllmaras< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Reimplemented in GRINS::SpalartAllmarasStabilizationBase< Viscosity >.

Definition at line 124 of file spalart_allmaras.C.

Referenced by GRINS::SpalartAllmarasStabilizationBase< Viscosity >::init_context().

125  {
126  // We should prerequest all the data
127  // we will need to build the linear system
128  // or evaluate a quantity of interest.
129  context.get_element_fe(_turbulence_vars.nu())->get_JxW();
130  context.get_element_fe(_turbulence_vars.nu())->get_phi();
131  context.get_element_fe(_turbulence_vars.nu())->get_dphi();
132  context.get_element_fe(_turbulence_vars.nu())->get_xyz();
133 
134  context.get_element_fe(_turbulence_vars.nu())->get_phi();
135  context.get_element_fe(_turbulence_vars.nu())->get_xyz();
136 
137  context.get_side_fe(_turbulence_vars.nu())->get_JxW();
138  context.get_side_fe(_turbulence_vars.nu())->get_phi();
139  context.get_side_fe(_turbulence_vars.nu())->get_dphi();
140  context.get_side_fe(_turbulence_vars.nu())->get_xyz();
141 
142  return;
143  }
TurbulenceFEVariables _turbulence_vars
template<class Mu >
void GRINS::SpalartAllmaras< Mu >::init_variables ( libMesh::FEMSystem *  system)
virtual

Initialize variables for this physics.

Reimplemented from GRINS::TurbulenceModelsBase< Viscosity >.

Reimplemented in GRINS::SpalartAllmarasStabilizationBase< Viscosity >, and GRINS::SpalartAllmarasSPGSMStabilization< Viscosity >.

Definition at line 94 of file spalart_allmaras.C.

References GRINS::TurbulenceModelsBase< Viscosity >::init_variables().

Referenced by GRINS::SpalartAllmarasStabilizationBase< Viscosity >::init_variables().

95  {
96  // Init base class.
98 
99  this->_turbulence_vars.init(system);
100  this->_flow_vars.init(system);
101  this->_press_var.init(system);
102 
103  // Init the variables belonging to SA helper
105 
106  // Initialize Boundary Mesh
107  this->boundary_mesh.reset(new libMesh::SerialMesh(system->get_mesh().comm() , this->_dim));
108 
109  // Use the _wall_ids set to build the boundary mesh object
110  (system->get_mesh()).boundary_info->sync(_wall_ids, *boundary_mesh);
111 
112  //this->distance_function.reset(new DistanceFunction(system->get_equation_systems(), dynamic_cast<libMesh::UnstructuredMesh&>(system->get_mesh()) ));
113  this->distance_function.reset(new DistanceFunction(system->get_equation_systems(), *boundary_mesh));
114 
115  // For now, we are hacking this. Without this initialize function being called
116  // the distance variable will just be zero. For the channel flow, we are just
117  // going to analytically compute the wall distance
118  //this->distance_function->initialize();
119 
120  return;
121  }
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
libMesh::UniquePtr< DistanceFunction > distance_function
unsigned int _dim
Physical dimension of problem.
virtual void init(libMesh::FEMSystem *system)
Add variables to the system.
PressureFEVariable _press_var
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.
SpalartAllmarasHelper _spalart_allmaras_helper
TurbulenceFEVariables _turbulence_vars
std::set< libMesh::boundary_id_type > _wall_ids
VelocityFEVariables _flow_vars
libMesh::UniquePtr< libMesh::SerialMesh > boundary_mesh
void init_variables(libMesh::FEMSystem *system)
template<class K >
void GRINS::SpalartAllmaras< K >::mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part.

Reimplemented from GRINS::Physics.

Reimplemented in GRINS::SpalartAllmarasSPGSMStabilization< Viscosity >.

Definition at line 308 of file spalart_allmaras.C.

311  {
312 #ifdef GRINS_USE_GRVY_TIMERS
313  this->_timer->BeginTimer("SpalartAllmaras::mass_residual");
314 #endif
315 
316  // First we get some references to cell-specific data that
317  // will be used to assemble the linear system.
318 
319  // Element Jacobian * quadrature weights for interior integration.
320  const std::vector<libMesh::Real> &JxW =
321  context.get_element_fe(this->_turbulence_vars.nu())->get_JxW();
322 
323  // The viscosity shape functions at interior quadrature points.
324  const std::vector<std::vector<libMesh::Real> >& nu_phi =
325  context.get_element_fe(this->_turbulence_vars.nu())->get_phi();
326 
327  // The number of local degrees of freedom in each variable.
328  const unsigned int n_nu_dofs = context.get_dof_indices(this->_turbulence_vars.nu()).size();
329 
330  // The subvectors and submatrices we need to fill:
331  libMesh::DenseSubVector<libMesh::Real> &F = context.get_elem_residual(this->_turbulence_vars.nu());
332 
333  //libMesh::DenseSubMatrix<libMesh::Real> &M = context.get_elem_jacobian(this->_turbulence_vars.nu(), this->_turbulence_vars.nu());
334 
335  unsigned int n_qpoints = context.get_element_qrule().n_points();
336 
337  for (unsigned int qp = 0; qp != n_qpoints; ++qp)
338  {
339  // For the mass residual, we need to be a little careful.
340  // The time integrator is handling the time-discretization
341  // for us so we need to supply M(u_fixed)*u' for the residual.
342  // u_fixed will be given by the fixed_interior_value function
343  // while u' will be given by the interior_rate function.
344  libMesh::Real nu_dot;
345  context.interior_rate(this->_turbulence_vars.nu(), qp, nu_dot);
346 
347  for (unsigned int i = 0; i != n_nu_dofs; ++i)
348  {
349  F(i) += -JxW[qp]*this->_rho*nu_dot*nu_phi[i][qp];
350 
351  if( compute_jacobian )
352  {
353  libmesh_not_implemented();
354  }// End of check on Jacobian
355 
356  } // End of element dof loop
357 
358  } // End of the quadrature point loop
359 
360 #ifdef GRINS_USE_GRVY_TIMERS
361  this->_timer->EndTimer("SpalartAllmaras::mass_residual");
362 #endif
363 
364  return;
365  }
TurbulenceFEVariables _turbulence_vars
libMesh::Number _rho
Material parameters, read from input.
template<class Mu >
void GRINS::SpalartAllmaras< Mu >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiAccessor< libMesh::Number > &  param_pointer 
) const
virtual

Each subclass will register its copy of an independent.

Reimplemented from GRINS::TurbulenceModelsBase< Viscosity >.

Reimplemented in GRINS::SpalartAllmarasSPGSMStabilization< Viscosity >, and GRINS::SpalartAllmarasStabilizationBase< Viscosity >.

Definition at line 369 of file spalart_allmaras.C.

References GRINS::ParameterUser::register_parameter().

Referenced by GRINS::SpalartAllmarasStabilizationBase< Viscosity >::register_parameter().

372  {
373  ParameterUser::register_parameter(param_name, param_pointer);
374  this->_mu.register_parameter(param_name, param_pointer);
375  this->_sa_params.register_parameter(param_name, param_pointer);
376  }
SpalartAllmarasParameters _sa_params
Object handling the plethora of parameters.
Viscosity _mu
Viscosity object.
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiAccessor< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
template<class Mu >
void GRINS::SpalartAllmaras< Mu >::register_variables ( )
private

Definition at line 83 of file spalart_allmaras.C.

References GRINS::GRINSPrivate::VariableWarehouse::check_and_register_variable(), GRINS::VariablesParsing::pressure_section(), GRINS::VariablesParsing::turbulence_section(), and GRINS::VariablesParsing::velocity_section().

Referenced by GRINS::SpalartAllmaras< Viscosity >::SpalartAllmaras().

84  {
86  this->_press_var);
88  this->_flow_vars);
90  this->_turbulence_vars);
91  }
PressureFEVariable _press_var
static void check_and_register_variable(const std::string &var_name, const FEVariablesBase &variable)
First check if var_name is registered and then register.
static std::string velocity_section()
static std::string pressure_section()
TurbulenceFEVariables _turbulence_vars
VelocityFEVariables _flow_vars
static std::string turbulence_section()
template<class Mu >
void GRINS::SpalartAllmaras< Mu >::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtual

Sets velocity variables to be time-evolving.

Reimplemented from GRINS::Physics.

Definition at line 146 of file spalart_allmaras.C.

147  {
148  // Tell the system to march velocity forward in time, but
149  // leave p as a constraint only
150  system->time_evolving(this->_turbulence_vars.nu());
151 
152  return;
153  }
TurbulenceFEVariables _turbulence_vars

Member Data Documentation

template<class Viscosity >
VelocityFEVariables GRINS::SpalartAllmaras< Viscosity >::_flow_vars
protected

Definition at line 98 of file spalart_allmaras.h.

template<class Viscosity >
bool GRINS::SpalartAllmaras< Viscosity >::_infinite_distance
protected

Definition at line 116 of file spalart_allmaras.h.

template<class Viscosity >
unsigned int GRINS::SpalartAllmaras< Viscosity >::_no_of_walls
protected
template<class Viscosity >
PressureFEVariable GRINS::SpalartAllmaras< Viscosity >::_press_var
protected

Definition at line 99 of file spalart_allmaras.h.

template<class Viscosity >
SpalartAllmarasParameters GRINS::SpalartAllmaras< Viscosity >::_sa_params
protected

Object handling the plethora of parameters.

Definition at line 107 of file spalart_allmaras.h.

template<class Viscosity >
SpalartAllmarasHelper GRINS::SpalartAllmaras< Viscosity >::_spalart_allmaras_helper
protected

Definition at line 104 of file spalart_allmaras.h.

template<class Viscosity >
TurbulenceFEVariables GRINS::SpalartAllmaras< Viscosity >::_turbulence_vars
protected

Definition at line 101 of file spalart_allmaras.h.

template<class Viscosity >
std::set<libMesh::boundary_id_type> GRINS::SpalartAllmaras< Viscosity >::_wall_ids
protected
template<class Viscosity >
libMesh::UniquePtr<libMesh::SerialMesh> GRINS::SpalartAllmaras< Viscosity >::boundary_mesh

Definition at line 86 of file spalart_allmaras.h.

template<class Viscosity >
libMesh::UniquePtr<DistanceFunction> GRINS::SpalartAllmaras< Viscosity >::distance_function

Definition at line 83 of file spalart_allmaras.h.


The documentation for this class was generated from the following files:

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