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

Physics class for Velocity Penalty. More...

#include <velocity_penalty.h>

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

Public Member Functions

 VelocityPenalty (const std::string &physics_name, const GetPot &input)
 
 ~VelocityPenalty ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register postprocessing variables for VelocityPenalty. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
virtual void read_input_options (const GetPot &input)
 Read options from GetPot input file. More...
 
bool compute_force (const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialization of Navier-Stokes variables. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Sets velocity variables to be time-evolving. More...
 
virtual void register_parameter (const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
 Each subclass will register its copy of an independent. More...
 
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 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 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 nonlocal_mass_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Mass matrix part(s) for scalar variables. More...
 
void init_bcs (libMesh::FEMSystem *system)
 
void init_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
void attach_neumann_bound_func (GRINS::NBCContainer &neumann_bcs)
 
void attach_dirichlet_bound_func (const GRINS::DBCContainer &dirichlet_bc)
 
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_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
BCHandlingBaseget_bc_handler ()
 
ICHandlingBaseget_ic_handler ()
 
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...
 

Protected Attributes

bool _quadratic_scaling
 
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > normal_vector_function
 
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > base_velocity_function
 
unsigned int _dim
 Physical dimension of problem. More...
 
PrimitiveFlowFEVariables _flow_vars
 
libMesh::Number _rho
 Material parameters, read from input. More...
 
Viscosity _mu
 Viscosity object. More...
 
const PhysicsName _physics_name
 Name of the physics object. Used for reading physics specific inputs. More...
 
GRINS::BCHandlingBase_bc_handler
 
GRINS::ICHandlingBase_ic_handler
 
std::set< libMesh::subdomain_id_type > _enabled_subdomains
 Subdomains on which the current Physics class is enabled. More...
 
bool _is_axisymmetric
 

Static Protected Attributes

static bool _is_steady = false
 Caches whether or not the solver that's being used is steady or not. More...
 

Private Member Functions

 VelocityPenalty ()
 

Private Attributes

unsigned int _velocity_penalty_x_index
 Index from registering this quantity. More...
 
unsigned int _velocity_penalty_y_index
 Index from registering this quantity. More...
 
unsigned int _velocity_penalty_z_index
 Index from registering this quantity. More...
 
unsigned int _velocity_penalty_base_x_index
 Index from registering this quantity. More...
 
unsigned int _velocity_penalty_base_y_index
 Index from registering this quantity. More...
 
unsigned int _velocity_penalty_base_z_index
 Index from registering this quantity. More...
 

Detailed Description

template<class Viscosity>
class GRINS::VelocityPenalty< Viscosity >

Physics class for Velocity Penalty.

Definition at line 50 of file velocity_penalty.h.

Constructor & Destructor Documentation

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

Definition at line 41 of file velocity_penalty.C.

42  : VelocityPenaltyBase<Mu>(physics_name, input),
49  {
50  return;
51  }
unsigned int _velocity_penalty_base_y_index
Index from registering this quantity.
unsigned int _velocity_penalty_y_index
Index from registering this quantity.
unsigned int _velocity_penalty_x_index
Index from registering this quantity.
unsigned int _velocity_penalty_base_x_index
Index from registering this quantity.
unsigned int _velocity_penalty_base_z_index
Index from registering this quantity.
unsigned int _velocity_penalty_z_index
Index from registering this quantity.
template<class Mu >
GRINS::VelocityPenalty< Mu >::~VelocityPenalty ( )

Definition at line 54 of file velocity_penalty.C.

55  {
56  return;
57  }
template<class Viscosity >
GRINS::VelocityPenalty< Viscosity >::VelocityPenalty ( )
private

Member Function Documentation

void GRINS::Physics::attach_dirichlet_bound_func ( const GRINS::DBCContainer dirichlet_bc)
inherited

Definition at line 150 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_dirichlet_bound_func().

151  {
152  _bc_handler->attach_dirichlet_bound_func( dirichlet_bc );
153  return;
154  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_dirichlet_bound_func(const GRINS::DBCContainer &dirichlet_bc)
void GRINS::Physics::attach_neumann_bound_func ( GRINS::NBCContainer neumann_bcs)
inherited

Definition at line 144 of file physics.C.

References GRINS::Physics::_bc_handler, and GRINS::BCHandlingBase::attach_neumann_bound_func().

145  {
146  _bc_handler->attach_neumann_bound_func( neumann_bcs );
147  return;
148  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
void attach_neumann_bound_func(GRINS::NBCContainer &neumann_bcs)
void GRINS::Physics::auxiliary_init ( MultiphysicsSystem system)
virtualinherited

Any auxillary initialization a Physics class may need.

This is called after all variables are added, so this method can safely query the MultiphysicsSystem about variable information.

Definition at line 113 of file physics.C.

114  {
115  return;
116  }
void GRINS::Physics::compute_element_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 185 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::element_constraint().

187  {
188  return;
189  }
void GRINS::Physics::compute_element_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited
template<class Mu >
bool GRINS::VelocityPenaltyBase< Mu >::compute_force ( const libMesh::Point &  point,
const libMesh::Real  time,
const libMesh::NumberVectorValue &  U,
libMesh::NumberVectorValue &  F,
libMesh::NumberTensorValue *  dFdU = NULL 
)
inherited

Definition at line 100 of file velocity_penalty_base.C.

105  {
106  // Velocity discrepancy (current velocity minus base velocity)
107  // normal to constraint plane, scaled by constraint penalty
108  // value
109  libmesh_assert(normal_vector_function.get());
110  libmesh_assert(base_velocity_function.get());
111 
112  libMesh::DenseVector<libMesh::Number> output_vec(3);
113 
114  (*normal_vector_function)(point, time,
115  output_vec);
116 
117  libMesh::NumberVectorValue U_N(output_vec(0),
118  output_vec(1),
119  output_vec(2));
120 
121  (*base_velocity_function)(point, time,
122  output_vec);
123 
124  const libMesh::NumberVectorValue U_B(output_vec(0),
125  output_vec(1),
126  output_vec(2));
127 
128  const libMesh::NumberVectorValue U_Rel = U-U_B;
129 
130  // Old code
131  // const libMesh::NumberVectorValue F1 = (U_Rel*U_N)*U_N; //
132 
133  // With correct sign and more natural normalization
134  const libMesh::Number U_N_mag = std::sqrt(U_N*U_N);
135 
136  if (!U_N_mag)
137  return false;
138 
139  const libMesh::NumberVectorValue U_N_unit = U_N/U_N_mag;
140 
141  F = -(U_Rel*U_N)*U_N_unit;
142 
143  if (dFdU)
144  for (unsigned int i=0; i != 3; ++i)
145  for (unsigned int j=0; j != 3; ++j)
146  (*dFdU)(i,j) = -(U_N(j))*U_N_unit(i);
147 
148  // With quadratic scaling
149  if (_quadratic_scaling)
150  {
151  const libMesh::Number U_Rel_mag = std::sqrt(U_Rel * U_Rel);
152 
153  // Modify dFdU first so as to reuse the old value of F
154  if (dFdU)
155  {
156  // dU_Rel/dU = I
157  // d(U_Rel*U_Rel)/dU = 2*U_Rel
158  // d|U_Rel|/dU = U_Rel/|U_Rel|
159 
160  const libMesh::NumberVectorValue U_Rel_unit = U_Rel/U_Rel_mag;
161 
162  (*dFdU) *= U_Rel_mag;
163 
164  if (U_Rel_mag)
165  for (unsigned int i=0; i != 3; ++i)
166  for (unsigned int j=0; j != 3; ++j)
167  (*dFdU)(i,j) += F(i)*U_Rel_unit(j);
168  }
169 
170  F *= U_Rel_mag;
171  }
172 
173  // With correction term to avoid doing work on flow
174  bool do_correction = false;
175  if (do_correction)
176  {
177  if (dFdU)
178  libmesh_not_implemented();
179 
180  const libMesh::Number U_Rel_mag_sq = U_Rel * U_Rel;
181  if (U_Rel_mag_sq)
182  {
183  F -= (U_Rel*F)*U_Rel/U_Rel_mag_sq;
184  }
185  }
186  return true;
187  }
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > normal_vector_function
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > base_velocity_function
void GRINS::Physics::compute_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 203 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::mass_residual().

205  {
206  return;
207  }
void GRINS::Physics::compute_nonlocal_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 197 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

199  {
200  return;
201  }
void GRINS::Physics::compute_nonlocal_mass_residual_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 209 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

211  {
212  return;
213  }
void GRINS::Physics::compute_nonlocal_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 179 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

181  {
182  return;
183  }
template<class Mu >
void GRINS::VelocityPenalty< Mu >::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtual

Reimplemented from GRINS::Physics.

Definition at line 242 of file velocity_penalty.C.

246  {
247  libMesh::DenseVector<libMesh::Number> output_vec(3);
248 
249  if( quantity_index == this->_velocity_penalty_x_index )
250  {
251  (*this->normal_vector_function)(point, context.time, output_vec);
252 
253  value = output_vec(0);
254  }
255  else if( quantity_index == this->_velocity_penalty_y_index )
256  {
257  (*this->normal_vector_function)(point, context.time, output_vec);
258 
259  value = output_vec(1);
260  }
261  else if( quantity_index == this->_velocity_penalty_z_index )
262  {
263  (*this->normal_vector_function)(point, context.time, output_vec);
264 
265  value = output_vec(2);
266  }
267  else if( quantity_index == this->_velocity_penalty_base_x_index )
268  {
269  (*this->base_velocity_function)(point, context.time, output_vec);
270 
271  value = output_vec(0);
272  }
273  else if( quantity_index == this->_velocity_penalty_base_y_index )
274  {
275  (*this->base_velocity_function)(point, context.time, output_vec);
276 
277  value = output_vec(1);
278  }
279  else if( quantity_index == this->_velocity_penalty_base_z_index )
280  {
281  (*this->base_velocity_function)(point, context.time, output_vec);
282 
283  value = output_vec(2);
284  }
285 
286  return;
287  }
unsigned int _velocity_penalty_base_y_index
Index from registering this quantity.
unsigned int _velocity_penalty_y_index
Index from registering this quantity.
unsigned int _velocity_penalty_x_index
Index from registering this quantity.
unsigned int _velocity_penalty_base_x_index
Index from registering this quantity.
unsigned int _velocity_penalty_base_z_index
Index from registering this quantity.
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > normal_vector_function
unsigned int _velocity_penalty_z_index
Index from registering this quantity.
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > base_velocity_function
void GRINS::Physics::compute_side_constraint_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Definition at line 191 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

193  {
194  return;
195  }
void GRINS::Physics::compute_side_time_derivative_cache ( const AssemblyContext context,
CachedValues cache 
)
virtualinherited

Reimplemented in GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 173 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_time_derivative().

175  {
176  return;
177  }
void GRINS::Physics::element_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited
template<class Mu >
void GRINS::VelocityPenalty< 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.

Definition at line 126 of file velocity_penalty.C.

129  {
130 #ifdef GRINS_USE_GRVY_TIMERS
131  this->_timer->BeginTimer("VelocityPenalty::element_time_derivative");
132 #endif
133 
134  // Element Jacobian * quadrature weights for interior integration
135  const std::vector<libMesh::Real> &JxW =
136  context.get_element_fe(this->_flow_vars.u_var())->get_JxW();
137 
138  // The shape functions at interior quadrature points.
139  const std::vector<std::vector<libMesh::Real> >& u_phi =
140  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
141 
142  const std::vector<libMesh::Point>& u_qpoint =
143  context.get_element_fe(this->_flow_vars.u_var())->get_xyz();
144 
145  // The number of local degrees of freedom in each variable
146  const unsigned int n_u_dofs = context.get_dof_indices(this->_flow_vars.u_var()).size();
147 
148  // The subvectors and submatrices we need to fill:
149  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.u_var()); // R_{u},{u}
150  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.v_var()); // R_{u},{v}
151  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.u_var()); // R_{v},{u}
152  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.v_var()); // R_{v},{v}
153 
154  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
155  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
156  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
157  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
158  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
159 
160  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_flow_vars.u_var()); // R_{u}
161  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_flow_vars.v_var()); // R_{v}
162  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
163 
164  if( this->_dim == 3 )
165  {
166  Kuw = &context.get_elem_jacobian(this->_flow_vars.u_var(), this->_flow_vars.w_var()); // R_{u},{w}
167  Kvw = &context.get_elem_jacobian(this->_flow_vars.v_var(), this->_flow_vars.w_var()); // R_{v},{w}
168 
169  Kwu = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.u_var()); // R_{w},{u}
170  Kwv = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.v_var()); // R_{w},{v}
171  Kww = &context.get_elem_jacobian(this->_flow_vars.w_var(), this->_flow_vars.w_var()); // R_{w},{w}
172  Fw = &context.get_elem_residual(this->_flow_vars.w_var()); // R_{w}
173  }
174 
175  unsigned int n_qpoints = context.get_element_qrule().n_points();
176 
177  for (unsigned int qp=0; qp != n_qpoints; qp++)
178  {
179  // Compute the solution at the old Newton iterate.
180  libMesh::Number u, v;
181  u = context.interior_value(this->_flow_vars.u_var(), qp);
182  v = context.interior_value(this->_flow_vars.v_var(), qp);
183 
184  libMesh::NumberVectorValue U(u,v);
185  if (this->_dim == 3)
186  U(2) = context.interior_value(this->_flow_vars.w_var(), qp); // w
187 
188  libMesh::NumberVectorValue F;
189  libMesh::NumberTensorValue dFdU;
190  libMesh::NumberTensorValue* dFdU_ptr =
191  compute_jacobian ? &dFdU : NULL;
192  if (!this->compute_force(u_qpoint[qp], context.time, U, F, dFdU_ptr))
193  continue;
194 
195  const libMesh::Real jac = JxW[qp];
196 
197  for (unsigned int i=0; i != n_u_dofs; i++)
198  {
199  const libMesh::Number jac_i = jac * u_phi[i][qp];
200 
201  Fu(i) += F(0)*jac_i;
202 
203  Fv(i) += F(1)*jac_i;
204  if( this->_dim == 3 )
205  {
206  (*Fw)(i) += F(2)*jac_i;
207  }
208 
209  if( compute_jacobian )
210  {
211  for (unsigned int j=0; j != n_u_dofs; j++)
212  {
213  const libMesh::Number jac_ij = context.get_elem_solution_derivative() * jac_i * u_phi[j][qp];
214  Kuu(i,j) += jac_ij * dFdU(0,0);
215  Kuv(i,j) += jac_ij * dFdU(0,1);
216  Kvu(i,j) += jac_ij * dFdU(1,0);
217  Kvv(i,j) += jac_ij * dFdU(1,1);
218 
219  if( this->_dim == 3 )
220  {
221  (*Kuw)(i,j) += jac_ij * dFdU(0,2);
222  (*Kvw)(i,j) += jac_ij * dFdU(1,2);
223 
224  (*Kwu)(i,j) += jac_ij * dFdU(2,0);
225  (*Kwv)(i,j) += jac_ij * dFdU(2,1);
226  (*Kww)(i,j) += jac_ij * dFdU(2,2);
227  }
228  }
229  }
230  }
231  }
232 
233 
234 #ifdef GRINS_USE_GRVY_TIMERS
235  this->_timer->EndTimer("VelocityPenalty::element_time_derivative");
236 #endif
237 
238  return;
239  }
unsigned int _dim
Physical dimension of problem.
bool compute_force(const libMesh::Point &point, const libMesh::Real time, const libMesh::NumberVectorValue &U, libMesh::NumberVectorValue &F, libMesh::NumberTensorValue *dFdU=NULL)
bool GRINS::Physics::enabled_on_elem ( const libMesh::Elem *  elem)
virtualinherited

Find if current physics is active on supplied element.

Definition at line 83 of file physics.C.

References GRINS::Physics::_enabled_subdomains.

84  {
85  // Check if enabled_subdomains flag has been set and if we're
86  // looking at a real element (rather than a nonlocal evaluation)
87  if( !elem || _enabled_subdomains.empty() )
88  return true;
89 
90  // Check if current physics is enabled on elem
91  if( _enabled_subdomains.find( elem->subdomain_id() ) == _enabled_subdomains.end() )
92  return false;
93 
94  return true;
95  }
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
BCHandlingBase * GRINS::Physics::get_bc_handler ( )
inlineinherited

Definition at line 282 of file physics.h.

References GRINS::Physics::_bc_handler.

283  {
284  return _bc_handler;
285  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
ICHandlingBase * GRINS::Physics::get_ic_handler ( )
inlineinherited

Definition at line 288 of file physics.h.

References GRINS::Physics::_ic_handler.

289  {
290  return _ic_handler;
291  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
void GRINS::Physics::init_bcs ( libMesh::FEMSystem *  system)
inherited

Definition at line 118 of file physics.C.

References GRINS::Physics::_bc_handler, GRINS::BCHandlingBase::init_bc_data(), GRINS::BCHandlingBase::init_dirichlet_bc_func_objs(), GRINS::BCHandlingBase::init_dirichlet_bcs(), and GRINS::BCHandlingBase::init_periodic_bcs().

119  {
120  // Only need to init BC's if the physics actually created a handler
121  if( _bc_handler )
122  {
123  _bc_handler->init_bc_data( *system );
124  _bc_handler->init_dirichlet_bcs( system );
126  _bc_handler->init_periodic_bcs( system );
127  }
128 
129  return;
130  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
virtual void init_dirichlet_bcs(libMesh::FEMSystem *system) const
virtual void init_periodic_bcs(libMesh::FEMSystem *system) const
virtual void init_dirichlet_bc_func_objs(libMesh::FEMSystem *system) const
virtual void init_bc_data(const libMesh::FEMSystem &system)
Override this method to initialize any system-dependent data.
template<class Mu >
void GRINS::VelocityPenalty< Mu >::init_context ( AssemblyContext context)
virtual

Initialize context for added physics variables.

Reimplemented from GRINS::IncompressibleNavierStokesBase< Viscosity >.

Definition at line 60 of file velocity_penalty.C.

61  {
62  context.get_element_fe(this->_flow_vars.u_var())->get_xyz();
63  context.get_element_fe(this->_flow_vars.u_var())->get_phi();
64 
65  return;
66  }
void GRINS::Physics::init_ics ( libMesh::FEMSystem *  system,
libMesh::CompositeFunction< libMesh::Number > &  all_ics 
)
inherited

Definition at line 133 of file physics.C.

References GRINS::Physics::_ic_handler, and GRINS::ICHandlingBase::init_ic_data().

135  {
136  if( _ic_handler )
137  {
138  _ic_handler->init_ic_data( *system, all_ics );
139  }
140 
141  return;
142  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
virtual void init_ic_data(const libMesh::FEMSystem &system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
Override this method to initialize any system-dependent data.
template<class Mu >
void GRINS::IncompressibleNavierStokesBase< Mu >::init_variables ( libMesh::FEMSystem *  system)
virtualinherited

Initialization of Navier-Stokes variables.

Add velocity and pressure variables to system.

Implements GRINS::Physics.

Reimplemented in GRINS::AveragedTurbineBase< Viscosity >, and GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >.

Definition at line 63 of file inc_navier_stokes_base.C.

Referenced by GRINS::IncompressibleNavierStokesStabilizationBase< Viscosity >::init_variables(), and GRINS::AveragedTurbineBase< Viscosity >::init_variables().

64  {
65  this->_dim = system->get_mesh().mesh_dimension();
66 
67  this->_flow_vars.init(system);
68 
69  this->_mu.init(system);
70 
71  return;
72  }
virtual void init(libMesh::FEMSystem *system)
unsigned int _dim
Physical dimension of problem.
bool GRINS::Physics::is_steady ( ) const
inherited

Returns whether or not this physics is being solved with a steady solver.

Definition at line 103 of file physics.C.

References GRINS::Physics::_is_steady.

Referenced by GRINS::Physics::set_is_steady().

104  {
105  return _is_steady;
106  }
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
void GRINS::Physics::mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited
void GRINS::Physics::nonlocal_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for scalar variables.

Reimplemented in GRINS::ScalarODE.

Definition at line 250 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_constraint().

253  {
254  return;
255  }
void GRINS::Physics::nonlocal_mass_residual ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Mass matrix part(s) for scalar variables.

Reimplemented in GRINS::ScalarODE, and GRINS::AveragedTurbine< Viscosity >.

Definition at line 264 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_mass_residual().

267  {
268  return;
269  }
void GRINS::Physics::nonlocal_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Time dependent part(s) of physics for scalar variables.

Reimplemented in GRINS::AveragedTurbine< Viscosity >, and GRINS::ScalarODE.

Definition at line 229 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::nonlocal_time_derivative().

232  {
233  return;
234  }
template<class Mu >
void GRINS::VelocityPenaltyBase< Mu >::read_input_options ( const GetPot &  input)
virtualinherited

Read options from GetPot input file.

Reimplemented from GRINS::Physics.

Definition at line 58 of file velocity_penalty_base.C.

References GRINS::velocity_penalty2, GRINS::velocity_penalty2_adjoint_stab, GRINS::velocity_penalty3, and GRINS::velocity_penalty3_adjoint_stab.

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

59  {
60  std::string base_physics_name = "VelocityPenalty";
61  if (this->_physics_name == velocity_penalty2 ||
63  base_physics_name.push_back('2');
64 
65  if (this->_physics_name == velocity_penalty3 ||
67  base_physics_name.push_back('3');
68 
69  std::string penalty_function =
70  input("Physics/"+base_physics_name+"/penalty_function",
71  std::string("0"));
72 
73  if (penalty_function == "0")
74  this->normal_vector_function.reset
75  (new libMesh::ZeroFunction<libMesh::Number>());
76  else
77  this->normal_vector_function.reset
78  (new libMesh::ParsedFunction<libMesh::Number>(penalty_function));
79 
80  std::string base_function =
81  input("Physics/"+base_physics_name+"/base_velocity",
82  std::string("0"));
83 
84  if (penalty_function == "0" && base_function == "0")
85  std::cout << "Warning! Zero VelocityPenalty specified!" << std::endl;
86 
87  if (base_function == "0")
88  this->base_velocity_function.reset
89  (new libMesh::ZeroFunction<libMesh::Number>());
90  else
91  this->base_velocity_function.reset
92  (new libMesh::ParsedFunction<libMesh::Number>(base_function));
93 
95  input("Physics/"+base_physics_name+"/quadratic_scaling", false);
96  }
const PhysicsName velocity_penalty2_adjoint_stab
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > normal_vector_function
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:254
const PhysicsName velocity_penalty3
const PhysicsName velocity_penalty2
libMesh::AutoPtr< libMesh::FunctionBase< libMesh::Number > > base_velocity_function
const PhysicsName velocity_penalty3_adjoint_stab
template<class Mu >
void GRINS::IncompressibleNavierStokesBase< Mu >::register_parameter ( const std::string &  param_name,
libMesh::ParameterMultiPointer< libMesh::Number > &  param_pointer 
) const
virtualinherited

Each subclass will register its copy of an independent.

Reimplemented from GRINS::ParameterUser.

Definition at line 114 of file inc_navier_stokes_base.C.

References GRINS::ParameterUser::register_parameter().

117  {
118  ParameterUser::register_parameter(param_name, param_pointer);
119  _mu.register_parameter(param_name, param_pointer);
120  }
virtual void register_parameter(const std::string &param_name, libMesh::ParameterMultiPointer< libMesh::Number > &param_pointer) const
Each subclass will register its copy of an independent.
template<class Mu >
void GRINS::VelocityPenalty< Mu >::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtual

Register postprocessing variables for VelocityPenalty.

Reimplemented from GRINS::Physics.

Definition at line 69 of file velocity_penalty.C.

References GRINS::PostProcessedQuantities< NumericType >::register_quantity().

71  {
72  std::string section = "Physics/"+this->_physics_name+"/output_vars";
73 
74  std::string vel_penalty = "vel_penalty";
75  if (this->_physics_name == "VelocityPenalty2")
76  vel_penalty += '2';
77 
78  if (this->_physics_name == "VelocityPenalty3")
79  vel_penalty += '3';
80 
81  if( input.have_variable(section) )
82  {
83  unsigned int n_vars = input.vector_variable_size(section);
84 
85  for( unsigned int v = 0; v < n_vars; v++ )
86  {
87  std::string name = input(section,"DIE!",v);
88 
89  if( name == std::string("velocity_penalty") )
90  {
92  postprocessing.register_quantity( vel_penalty+"_x" );
93 
95  postprocessing.register_quantity( vel_penalty+"_y" );
96 
98  postprocessing.register_quantity( vel_penalty+"_z" );
99  }
100  else if( name == std::string("velocity_penalty_base") )
101  {
103  postprocessing.register_quantity( vel_penalty+"_base_x" );
104 
106  postprocessing.register_quantity( vel_penalty+"_base_y" );
107 
109  postprocessing.register_quantity( vel_penalty+"_base_z" );
110  }
111  else
112  {
113  std::cerr << "Error: Invalue output_vars value for "+this->_physics_name << std::endl
114  << " Found " << name << std::endl
115  << " Acceptable values are: velocity_penalty" << std::endl
116  << " velocity_penalty_base" << std::endl;
117  libmesh_error();
118  }
119  }
120  }
121 
122  return;
123  }
unsigned int _velocity_penalty_base_y_index
Index from registering this quantity.
unsigned int _velocity_penalty_y_index
Index from registering this quantity.
unsigned int _velocity_penalty_x_index
Index from registering this quantity.
unsigned int _velocity_penalty_base_x_index
Index from registering this quantity.
unsigned int _velocity_penalty_base_z_index
Index from registering this quantity.
unsigned int _velocity_penalty_z_index
Index from registering this quantity.
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:254
void GRINS::Physics::set_is_steady ( bool  is_steady)
inherited

Sets whether this physics is to be solved with a steady solver or not.

Since the member variable is static, only needs to be called on a single physics.

Definition at line 97 of file physics.C.

References GRINS::Physics::_is_steady, and GRINS::Physics::is_steady().

98  {
100  return;
101  }
bool is_steady() const
Returns whether or not this physics is being solved with a steady solver.
Definition: physics.C:103
static bool _is_steady
Caches whether or not the solver that's being used is steady or not.
Definition: physics.h:266
void GRINS::ParameterUser::set_parameter ( libMesh::Number &  param_variable,
const GetPot &  input,
const std::string &  param_name,
libMesh::Number  param_default 
)
virtualinherited

Each subclass can simultaneously read a parameter value from.

Definition at line 35 of file parameter_user.C.

References GRINS::ParameterUser::_my_name, and GRINS::ParameterUser::_my_parameters.

Referenced by GRINS::AveragedFanAdjointStabilization< Viscosity >::AveragedFanAdjointStabilization(), GRINS::AveragedTurbineAdjointStabilization< Viscosity >::AveragedTurbineAdjointStabilization(), GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::BoussinesqBuoyancyAdjointStabilization(), GRINS::BoussinesqBuoyancyBase::BoussinesqBuoyancyBase(), GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::BoussinesqBuoyancySPGSMStabilization(), GRINS::ConstantConductivity::ConstantConductivity(), GRINS::ConstantPrandtlConductivity::ConstantPrandtlConductivity(), GRINS::ConstantSourceFunction::ConstantSourceFunction(), GRINS::ConstantSourceTerm::ConstantSourceTerm(), GRINS::ConstantSpecificHeat::ConstantSpecificHeat(), GRINS::ConstantViscosity::ConstantViscosity(), GRINS::ElasticCable< StressStrainLaw >::ElasticCable(), GRINS::ElasticCableConstantGravity::ElasticCableConstantGravity(), GRINS::ElasticMembrane< StressStrainLaw >::ElasticMembrane(), GRINS::ElasticMembraneConstantPressure::ElasticMembraneConstantPressure(), GRINS::HeatConduction< Conductivity >::HeatConduction(), GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::IncompressibleNavierStokesBase< Viscosity >::IncompressibleNavierStokesBase(), GRINS::AverageNusseltNumber::init(), GRINS::MooneyRivlin::MooneyRivlin(), GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >::ReactingLowMachNavierStokesBase(), GRINS::HookesLaw1D::read_input_options(), GRINS::HookesLaw::read_input_options(), GRINS::AxisymmetricBoussinesqBuoyancy::read_input_options(), and GRINS::VelocityDragAdjointStabilization< Viscosity >::VelocityDragAdjointStabilization().

39  {
40  param_variable = input(param_name, param_default);
41 
42  libmesh_assert_msg(!_my_parameters.count(param_name),
43  "ERROR: " << _my_name << " double-registered parameter " <<
44  param_name);
45 
46  _my_parameters[param_name] = &param_variable;
47  }
std::map< std::string, libMesh::Number * > _my_parameters
template<class Mu >
void GRINS::IncompressibleNavierStokesBase< Mu >::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtualinherited

Sets velocity variables to be time-evolving.

Reimplemented from GRINS::Physics.

Reimplemented in GRINS::AveragedTurbineBase< Viscosity >, and GRINS::ParsedVelocitySourceBase< Viscosity >.

Definition at line 75 of file inc_navier_stokes_base.C.

Referenced by GRINS::AveragedTurbineBase< Viscosity >::set_time_evolving_vars().

76  {
77  const unsigned int dim = system->get_mesh().mesh_dimension();
78 
79  // Tell the system to march velocity forward in time, but
80  // leave p as a constraint only
81  system->time_evolving(_flow_vars.u_var());
82  system->time_evolving(_flow_vars.v_var());
83 
84  if (dim == 3)
85  system->time_evolving(_flow_vars.w_var());
86 
87  return;
88  }
void GRINS::Physics::side_constraint ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Constraint part(s) of physics for boundaries of elements on the domain boundary.

Reimplemented in GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 243 of file physics.C.

Referenced by GRINS::MultiphysicsSystem::side_constraint().

246  {
247  return;
248  }
void GRINS::Physics::side_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtualinherited

Member Data Documentation

GRINS::BCHandlingBase* GRINS::Physics::_bc_handler
protectedinherited
template<class Viscosity >
unsigned int GRINS::IncompressibleNavierStokesBase< Viscosity >::_dim
protectedinherited

Physical dimension of problem.

Todo:
Do we really need to cache this?

Definition at line 79 of file inc_navier_stokes_base.h.

std::set<libMesh::subdomain_id_type> GRINS::Physics::_enabled_subdomains
protectedinherited

Subdomains on which the current Physics class is enabled.

Definition at line 261 of file physics.h.

Referenced by GRINS::Physics::enabled_on_elem(), and GRINS::Physics::read_input_options().

template<class Viscosity >
PrimitiveFlowFEVariables GRINS::IncompressibleNavierStokesBase< Viscosity >::_flow_vars
protectedinherited

Definition at line 81 of file inc_navier_stokes_base.h.

GRINS::ICHandlingBase* GRINS::Physics::_ic_handler
protectedinherited
bool GRINS::Physics::_is_axisymmetric
protectedinherited
bool GRINS::Physics::_is_steady = false
staticprotectedinherited

Caches whether or not the solver that's being used is steady or not.

This is need, for example, in flow stabilization as the tau terms change depending on whether the solver is steady or unsteady.

Definition at line 266 of file physics.h.

Referenced by GRINS::Physics::is_steady(), and GRINS::Physics::set_is_steady().

template<class Viscosity >
Viscosity GRINS::IncompressibleNavierStokesBase< Viscosity >::_mu
protectedinherited

Viscosity object.

Definition at line 88 of file inc_navier_stokes_base.h.

const PhysicsName GRINS::Physics::_physics_name
protectedinherited

Name of the physics object. Used for reading physics specific inputs.

We use a reference because the physics names are const global objects in GRINS namespace

No, we use a copy, because otherwise as soon as the memory in std::set<std::string> requested_physics gets overwritten we get in trouble.

Definition at line 254 of file physics.h.

Referenced by GRINS::SourceTermBase::parse_var_info(), and GRINS::Physics::read_input_options().

template<class Viscosity >
bool GRINS::VelocityPenaltyBase< Viscosity >::_quadratic_scaling
protectedinherited

Definition at line 64 of file velocity_penalty_base.h.

template<class Viscosity >
libMesh::Number GRINS::IncompressibleNavierStokesBase< Viscosity >::_rho
protectedinherited

Material parameters, read from input.

Todo:
Create objects to allow for function specification

Definition at line 85 of file inc_navier_stokes_base.h.

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

template<class Viscosity >
unsigned int GRINS::VelocityPenalty< Viscosity >::_velocity_penalty_base_x_index
private

Index from registering this quantity.

Definition at line 89 of file velocity_penalty.h.

template<class Viscosity >
unsigned int GRINS::VelocityPenalty< Viscosity >::_velocity_penalty_base_y_index
private

Index from registering this quantity.

Definition at line 92 of file velocity_penalty.h.

template<class Viscosity >
unsigned int GRINS::VelocityPenalty< Viscosity >::_velocity_penalty_base_z_index
private

Index from registering this quantity.

Definition at line 95 of file velocity_penalty.h.

template<class Viscosity >
unsigned int GRINS::VelocityPenalty< Viscosity >::_velocity_penalty_x_index
private

Index from registering this quantity.

Definition at line 80 of file velocity_penalty.h.

template<class Viscosity >
unsigned int GRINS::VelocityPenalty< Viscosity >::_velocity_penalty_y_index
private

Index from registering this quantity.

Definition at line 83 of file velocity_penalty.h.

template<class Viscosity >
unsigned int GRINS::VelocityPenalty< Viscosity >::_velocity_penalty_z_index
private

Index from registering this quantity.

Definition at line 86 of file velocity_penalty.h.

template<class Viscosity >
libMesh::AutoPtr<libMesh::FunctionBase<libMesh::Number> > GRINS::VelocityPenaltyBase< Viscosity >::base_velocity_function
protectedinherited

Definition at line 68 of file velocity_penalty_base.h.

template<class Viscosity >
libMesh::AutoPtr<libMesh::FunctionBase<libMesh::Number> > GRINS::VelocityPenaltyBase< Viscosity >::normal_vector_function
protectedinherited

Definition at line 66 of file velocity_penalty_base.h.


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

Generated on Mon Jun 22 2015 21:32:24 for GRINS-0.6.0 by  doxygen 1.8.9.1