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

#include <elastic_cable.h>

Inheritance diagram for GRINS::ElasticCable< StressStrainLaw >:
Inheritance graph
[legend]
Collaboration diagram for GRINS::ElasticCable< StressStrainLaw >:
Collaboration graph
[legend]

Public Member Functions

 ElasticCable (const PhysicsName &physics_name, const GetPot &input, bool lambda_sq_var)
 
virtual ~ElasticCable ()
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register postprocessing variables for ElasticCable. 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 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 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 compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 Compute the registered postprocessed quantities. More...
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialize variables for this physics. More...
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Set which variables are time evolving. More...
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
virtual void read_input_options (const GetPot &input)
 Read options from GetPot input file. By default, nothing is read. 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 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 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...
 
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...
 

Protected Member Functions

const libMesh::FEGenericBase< libMesh::Real > * get_fe (const AssemblyContext &context)
 

Protected Attributes

SolidMechanicsFEVariables _disp_vars
 
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

 ElasticCable ()
 
libMesh::AutoPtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe (const libMesh::Elem &elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
 
void compute_metric_tensors (unsigned int qp, const libMesh::FEBase &elem, const AssemblyContext &context, const libMesh::Gradient &grad_u, const libMesh::Gradient &grad_v, const libMesh::Gradient &grad_w, libMesh::TensorValue< libMesh::Real > &a_cov, libMesh::TensorValue< libMesh::Real > &a_contra, libMesh::TensorValue< libMesh::Real > &A_cov, libMesh::TensorValue< libMesh::Real > &A_contra, libMesh::Real &lambda_sq)
 

Private Attributes

StressStrainLaw _stress_strain_law
 
libMesh::Real _A
 
bool _is_compressible
 
std::vector< unsigned int > _stress_indices
 Index from registering this quantity. Each component will have it's own index. More...
 
std::vector< unsigned int > _strain_indices
 Index from registering this quantity. Each component will have it's own index. More...
 
std::vector< unsigned int > _force_indices
 Index from registering this quantity. Each component will have it's own index. More...
 

Detailed Description

template<typename StressStrainLaw>
class GRINS::ElasticCable< StressStrainLaw >

Definition at line 38 of file elastic_cable.h.

Constructor & Destructor Documentation

template<typename StressStrainLaw >
GRINS::ElasticCable< StressStrainLaw >::ElasticCable ( const PhysicsName physics_name,
const GetPot &  input,
bool  lambda_sq_var 
)

Definition at line 46 of file elastic_cable.C.

References GRINS::ElasticCable< StressStrainLaw >::_A, GRINS::Physics::_bc_handler, GRINS::Physics::_ic_handler, and GRINS::ParameterUser::set_parameter().

48  : ElasticCableBase(physics_name,input),
49  _stress_strain_law(input),
50  _A( 1.0 ),
51  _is_compressible(is_compressible)
52  {
53  // Force the user to set A
54  if( !input.have_variable("Physics/"+physics_name+"/A") )
55  {
56  std::cerr << "Error: Must specify initial area for "+physics_name << std::endl
57  << " Input the option Physics/"+physics_name+"/A" << std::endl;
58  libmesh_error();
59  }
60 
61  this->set_parameter
62  (_A, input, "Physics/"+physics_name+"/A", _A );
63 
64  this->_bc_handler = new SolidMechanicsBCHandling( physics_name, input );
65 
66  this->_ic_handler = new GenericICHandler(physics_name, input);
67 
68  return;
69  }
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.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
libMesh::Real _A
Definition: elastic_cable.h:93
StressStrainLaw _stress_strain_law
Definition: elastic_cable.h:91
template<typename StressStrainLaw >
GRINS::ElasticCable< StressStrainLaw >::~ElasticCable ( )
virtual

Definition at line 72 of file elastic_cable.C.

73  {
74  return;
75  }
template<typename StressStrainLaw >
GRINS::ElasticCable< StressStrainLaw >::ElasticCable ( )
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  }
template<typename StressStrainLaw >
libMesh::AutoPtr< libMesh::FEGenericBase< libMesh::Real > > GRINS::ElasticCable< StressStrainLaw >::build_new_fe ( const libMesh::Elem &  elem,
const libMesh::FEGenericBase< libMesh::Real > *  fe,
const libMesh::Point  p 
)
private

Definition at line 414 of file elastic_cable.C.

417  {
418  using namespace libMesh;
419  FEType fe_type = fe->get_fe_type();
420 
421  // If we don't have an Elem to evaluate on, then the only functions
422  // we can sensibly evaluate are the scalar dofs which are the same
423  // everywhere.
424  libmesh_assert(&elem || fe_type.family == SCALAR);
425 
426  unsigned int elem_dim = &elem ? elem.dim() : 0;
427 
428  AutoPtr<FEGenericBase<libMesh::Real> >
429  fe_new(FEGenericBase<libMesh::Real>::build(elem_dim, fe_type));
430 
431  // Map the physical co-ordinates to the master co-ordinates using the inverse_map from fe_interface.h
432  // Build a vector of point co-ordinates to send to reinit
433  Point master_point = &elem ?
434  FEInterface::inverse_map(elem_dim, fe_type, &elem, p) :
435  Point(0);
436 
437  std::vector<Point> coor(1, master_point);
438 
439  // Reinitialize the element and compute the shape function values at coor
440  fe_new->reinit (&elem, &coor);
441 
442  return fe_new;
443  }
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
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  }
template<typename StressStrainLaw >
void GRINS::ElasticCable< StressStrainLaw >::compute_metric_tensors ( unsigned int  qp,
const libMesh::FEBase &  elem,
const AssemblyContext context,
const libMesh::Gradient &  grad_u,
const libMesh::Gradient &  grad_v,
const libMesh::Gradient &  grad_w,
libMesh::TensorValue< libMesh::Real > &  a_cov,
libMesh::TensorValue< libMesh::Real > &  a_contra,
libMesh::TensorValue< libMesh::Real > &  A_cov,
libMesh::TensorValue< libMesh::Real > &  A_contra,
libMesh::Real &  lambda_sq 
)
private

Definition at line 446 of file elastic_cable.C.

457  {
458  const std::vector<libMesh::RealGradient>& dxdxi = elem.get_dxyzdxi();
459 
460  const std::vector<libMesh::Real>& dxidx = elem.get_dxidx();
461  const std::vector<libMesh::Real>& dxidy = elem.get_dxidy();
462  const std::vector<libMesh::Real>& dxidz = elem.get_dxidz();
463 
464  libMesh::RealGradient dxi( dxidx[qp], dxidy[qp], dxidz[qp] );
465 
466  libMesh::RealGradient dudxi( grad_u(0), grad_v(0), grad_w(0) );
467 
468  // Covariant metric tensor of reference configuration
469  a_cov.zero();
470  a_cov(0,0) = dxdxi[qp]*dxdxi[qp];
471  a_cov(1,1) = 1.0;
472  a_cov(2,2) = 1.0;
473 
474  // Covariant metric tensor of current configuration
475  A_cov.zero();
476  A_cov(0,0) = (dxdxi[qp] + dudxi)*(dxdxi[qp] + dudxi);
477 
478  // Contravariant metric tensor of reference configuration
479  a_contra.zero();
480  a_contra(0,0) = 1/a_cov(0,0);
481  a_contra(1,1) = 1.0;
482  a_contra(2,2) = 1.0;
483 
484  // Contravariant metric tensor in current configuration is A_cov^{-1}
485  A_contra.zero();
486  A_contra(0,0) = 1/A_cov(0,0);
487 
488  // If the material is compressible, then lambda_sq is an independent variable
489  if( _is_compressible )
490  {
491  libmesh_not_implemented();
492  //lambda_sq = context.interior_value(this->_lambda_sq_var, qp);
493  }
494  else
495  {
496  // If the material is incompressible, lambda^2 is known
497  lambda_sq = a_cov(0,0)/A_cov(0,0);//det_a/det_A;
498  }
499 
500  //Update the covariant and contravariant tensors of current configuration
501  A_cov(1,1) = lambda_sq;
502  A_cov(2,2) = lambda_sq;
503  A_contra(1,1) = 1.0/lambda_sq;
504  A_contra(2,2) = 1.0/lambda_sq;
505 
506  return;
507  }
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<typename StressStrainLaw >
void GRINS::ElasticCable< StressStrainLaw >::compute_postprocessed_quantity ( unsigned int  quantity_index,
const AssemblyContext context,
const libMesh::Point &  point,
libMesh::Real &  value 
)
virtual

Compute the registered postprocessed quantities.

Reimplemented from GRINS::Physics.

Definition at line 308 of file elastic_cable.C.

312  {
313  bool is_strain = false;
314  if( !_strain_indices.empty() )
315  is_strain = ( _strain_indices[0] == quantity_index );
316 
317  bool is_stress = false;
318  if( !_stress_indices.empty() )
319  is_stress = ( _stress_indices[0] == quantity_index );
320 
321  bool is_force = false;
322  if( !_force_indices.empty() )
323  is_force = ( _force_indices[0] == quantity_index );
324 
325  if( is_strain || is_stress || is_force )
326  {
327  const unsigned int n_u_dofs = context.get_dof_indices(_disp_vars.u_var()).size();
328 
329  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( _disp_vars.u_var() );
330  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( _disp_vars.v_var() );
331  const libMesh::DenseSubVector<libMesh::Number>& w_coeffs = context.get_elem_solution( _disp_vars.w_var() );
332 
333  // Build new FE for the current point. We need this to build tensors at point.
334  libMesh::AutoPtr<libMesh::FEGenericBase<libMesh::Real> > fe_new = this->build_new_fe( context.get_elem(), this->get_fe(context), point );
335 
336  const std::vector<std::vector<libMesh::Real> >& dphi_dxi = fe_new->get_dphidxi();
337 
338  // Need these to build up the covariant and contravariant metric tensors
339  const std::vector<libMesh::RealGradient>& dxdxi = fe_new->get_dxyzdxi();
340 
341  // Gradients are w.r.t. master element coordinates
342  libMesh::Gradient grad_u, grad_v, grad_w;
343  for( unsigned int d = 0; d < n_u_dofs; d++ )
344  {
345  libMesh::RealGradient u_gradphi( dphi_dxi[d][0] );
346  grad_u += u_coeffs(d)*u_gradphi;
347  grad_v += v_coeffs(d)*u_gradphi;
348  grad_w += w_coeffs(d)*u_gradphi;
349  }
350 
351  libMesh::RealGradient grad_x( dxdxi[0](0) );
352  libMesh::RealGradient grad_y( dxdxi[0](1) );
353  libMesh::RealGradient grad_z( dxdxi[0](2) );
354 
355  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
356  libMesh::Real lambda_sq = 0;
357 
358  this->compute_metric_tensors(0, *fe_new, context, grad_u, grad_v, grad_w, a_cov, a_contra, A_cov, A_contra, lambda_sq );
359 
360  libMesh::Real det_a = a_cov(0,0)*a_cov(1,1) - a_cov(0,1)*a_cov(1,0);
361  libMesh::Real det_A = A_cov(0,0)*A_cov(1,1) - A_cov(0,1)*A_cov(1,0);
362 
363  libMesh::Real I3 = lambda_sq*det_A/det_a;
364 
366  _stress_strain_law.compute_stress(1,a_contra,a_cov,A_contra,A_cov,tau);
367 
368  // We have everything we need for strain now, so check if we are computing strain
369  if( is_strain )
370  {
371  if( _strain_indices[0] == quantity_index )
372  {
373  value = 0.5*(A_cov(0,0) - a_cov(0,0));
374  }
375  else
376  {
377  //Wat?!
378  libmesh_error();
379  }
380  return;
381  }
382 
383  if( is_stress )
384  {
385  if( _stress_indices[0] == quantity_index )
386  {
387  // Need to convert to Cauchy stress
388  value = tau(0,0)/std::sqrt(I3);
389  }
390  else
391  {
392  libmesh_error();
393  }
394  }
395  if( is_force )
396  {
397 
398  if( _force_indices[0] == quantity_index )
399  {
400  //This force is in deformed configuration and will be significantly influenced by large strains
401  value = tau(0,0)/std::sqrt(I3)*_A;
402  }
403  else
404  {
405  libmesh_error();
406  }
407  }
408  }
409 
410  return;
411  }
void compute_metric_tensors(unsigned int qp, const libMesh::FEBase &elem, const AssemblyContext &context, const libMesh::Gradient &grad_u, const libMesh::Gradient &grad_v, const libMesh::Gradient &grad_w, libMesh::TensorValue< libMesh::Real > &a_cov, libMesh::TensorValue< libMesh::Real > &a_contra, libMesh::TensorValue< libMesh::Real > &A_cov, libMesh::TensorValue< libMesh::Real > &A_contra, libMesh::Real &lambda_sq)
void compute_stress(unsigned int dim, const libMesh::TensorValue< libMesh::Real > &g_contra, const libMesh::TensorValue< libMesh::Real > &g_cov, const libMesh::TensorValue< libMesh::Real > &G_contra, const libMesh::TensorValue< libMesh::Real > &G_cov, libMesh::TensorValue< libMesh::Real > &stress)
std::vector< unsigned int > _force_indices
Index from registering this quantity. Each component will have it's own index.
libMesh::Real _A
Definition: elastic_cable.h:93
SolidMechanicsFEVariables _disp_vars
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)
StressStrainLaw _stress_strain_law
Definition: elastic_cable.h:91
std::vector< unsigned int > _stress_indices
Index from registering this quantity. Each component will have it's own index.
Definition: elastic_cable.h:98
std::vector< unsigned int > _strain_indices
Index from registering this quantity. Each component will have it's own index.
libMesh::AutoPtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe(const libMesh::Elem &elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
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<typename StressStrainLaw >
void GRINS::ElasticCable< StressStrainLaw >::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 134 of file elastic_cable.C.

137  {
138  const unsigned int n_u_dofs = context.get_dof_indices(_disp_vars.u_var()).size();
139 
140  const std::vector<libMesh::Real> &JxW =
141  this->get_fe(context)->get_JxW();
142 
143  // Residuals that we're populating
144  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_disp_vars.u_var());
145  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_disp_vars.v_var());
146  libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(_disp_vars.w_var());
147 
148  //Grab the Jacobian matrix as submatrices
149  //libMesh::DenseMatrix<libMesh::Number> &K = context.get_elem_jacobian();
150  libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(_disp_vars.u_var(),_disp_vars.u_var());
151  libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(_disp_vars.u_var(),_disp_vars.v_var());
152  libMesh::DenseSubMatrix<libMesh::Number> &Kuw = context.get_elem_jacobian(_disp_vars.u_var(),_disp_vars.w_var());
153  libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(_disp_vars.v_var(),_disp_vars.u_var());
154  libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(_disp_vars.v_var(),_disp_vars.v_var());
155  libMesh::DenseSubMatrix<libMesh::Number> &Kvw = context.get_elem_jacobian(_disp_vars.v_var(),_disp_vars.w_var());
156  libMesh::DenseSubMatrix<libMesh::Number> &Kwu = context.get_elem_jacobian(_disp_vars.w_var(),_disp_vars.u_var());
157  libMesh::DenseSubMatrix<libMesh::Number> &Kwv = context.get_elem_jacobian(_disp_vars.w_var(),_disp_vars.v_var());
158  libMesh::DenseSubMatrix<libMesh::Number> &Kww = context.get_elem_jacobian(_disp_vars.w_var(),_disp_vars.w_var());
159 
160 
161  unsigned int n_qpoints = context.get_element_qrule().n_points();
162 
163  // All shape function gradients are w.r.t. master element coordinates
164  const std::vector<std::vector<libMesh::Real> >& dphi_dxi = this->get_fe(context)->get_dphidxi();
165 
166  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( _disp_vars.u_var() );
167  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( _disp_vars.v_var() );
168  const libMesh::DenseSubVector<libMesh::Number>& w_coeffs = context.get_elem_solution( _disp_vars.w_var() );
169 
170  // Need these to build up the covariant and contravariant metric tensors
171  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
172 
173  const unsigned int dim = 1; // The cable dimension is always 1 for this physics
174 
175 
176  for (unsigned int qp=0; qp != n_qpoints; qp++)
177  {
178  // Gradients are w.r.t. master element coordinates
179  libMesh::Gradient grad_u, grad_v, grad_w;
180 
181  for( unsigned int d = 0; d < n_u_dofs; d++ )
182  {
183  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp] );
184  grad_u += u_coeffs(d)*u_gradphi;
185  grad_v += v_coeffs(d)*u_gradphi;
186  grad_w += w_coeffs(d)*u_gradphi;
187  }
188 
189  libMesh::RealGradient grad_x( dxdxi[qp](0) );
190  libMesh::RealGradient grad_y( dxdxi[qp](1) );
191  libMesh::RealGradient grad_z( dxdxi[qp](2) );
192 
193  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
194  libMesh::Real lambda_sq = 0;
195 
196  this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
197  grad_u, grad_v, grad_w,
198  a_cov, a_contra, A_cov, A_contra,
199  lambda_sq );
200 
201  // Compute stress tensor
203  ElasticityTensor C;
204  //_stress_strain_law.compute_stress(dim,a_contra,a_cov,A_contra,A_cov,tau);
205  _stress_strain_law.compute_stress_and_elasticity(dim,a_contra,a_cov,A_contra,A_cov,tau,C);
206 
207 
208  libMesh::Real jac = JxW[qp];
209 
210  for (unsigned int i=0; i != n_u_dofs; i++)
211  {
212  libMesh::RealGradient u_gradphi( dphi_dxi[i][qp] );
213 
214  Fu(i) -= tau(0,0)*_A*( (grad_x(0) + grad_u(0))*u_gradphi(0) ) * jac;
215 
216  Fv(i) -= tau(0,0)*_A*( (grad_y(0) + grad_v(0))*u_gradphi(0) ) * jac;
217 
218  Fw(i) -= tau(0,0)*_A*( (grad_z(0) + grad_w(0))*u_gradphi(0) ) * jac;
219 
220  }
221 
222  if( compute_jacobian )
223  {
224  for(unsigned int i=0; i != n_u_dofs; i++)
225  {
226  libMesh::RealGradient u_gradphi_I( dphi_dxi[i][qp] );
227  for(unsigned int j=0; j != n_u_dofs; j++)
228  {
229  libMesh::RealGradient u_gradphi_J( dphi_dxi[j][qp] );
230 
231  const libMesh::Real diag_term = _A*jac*tau(0,0)*( u_gradphi_J(0)*u_gradphi_I(0));
232 
233  Kuu(i,j) -= diag_term;
234 
235  Kvv(i,j) -= diag_term;
236 
237  Kww(i,j) -= diag_term;
238 
239  const libMesh::Real dgamma_du = ( u_gradphi_J(0)*(grad_x(0)+grad_u(0)) );
240 
241  const libMesh::Real dgamma_dv = ( u_gradphi_J(0)*(grad_y(0)+grad_v(0)) );
242 
243  const libMesh::Real dgamma_dw = ( u_gradphi_J(0)*(grad_z(0)+grad_w(0)) );
244 
245  const libMesh::Real C1 = _A*jac*C(0,0,0,0);
246 
247  const libMesh::Real x_term = C1*( (grad_x(0)+grad_u(0))*u_gradphi_I(0) );
248 
249  const libMesh::Real y_term = C1*( (grad_y(0)+grad_v(0))*u_gradphi_I(0) );
250 
251  const libMesh::Real z_term = C1*( (grad_z(0)+grad_w(0))*u_gradphi_I(0) );
252 
253  Kuu(i,j) -= x_term*dgamma_du;
254 
255  Kuv(i,j) -= x_term*dgamma_dv;
256 
257  Kuw(i,j) -= x_term*dgamma_dw;
258 
259  Kvu(i,j) -= y_term*dgamma_du;
260 
261  Kvv(i,j) -= y_term*dgamma_dv;
262 
263  Kvw(i,j) -= y_term*dgamma_dw;
264 
265  Kwu(i,j) -= z_term*dgamma_du;
266 
267  Kwv(i,j) -= z_term*dgamma_dv;
268 
269  Kww(i,j) -= z_term*dgamma_dw;
270 
271  }
272  }
273  //libmesh_not_implemented();
274  }
275  }
276 
277  return;
278  }
void compute_metric_tensors(unsigned int qp, const libMesh::FEBase &elem, const AssemblyContext &context, const libMesh::Gradient &grad_u, const libMesh::Gradient &grad_v, const libMesh::Gradient &grad_w, libMesh::TensorValue< libMesh::Real > &a_cov, libMesh::TensorValue< libMesh::Real > &a_contra, libMesh::TensorValue< libMesh::Real > &A_cov, libMesh::TensorValue< libMesh::Real > &A_contra, libMesh::Real &lambda_sq)
libMesh::Real _A
Definition: elastic_cable.h:93
void compute_stress_and_elasticity(unsigned int dim, const libMesh::TensorValue< libMesh::Real > &g_contra, const libMesh::TensorValue< libMesh::Real > &g_cov, const libMesh::TensorValue< libMesh::Real > &G_contra, const libMesh::TensorValue< libMesh::Real > &G_cov, libMesh::TensorValue< libMesh::Real > &stress, ElasticityTensor &C)
SolidMechanicsFEVariables _disp_vars
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)
StressStrainLaw _stress_strain_law
Definition: elastic_cable.h:91
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
const libMesh::FEGenericBase< libMesh::Real > * GRINS::ElasticCableBase::get_fe ( const AssemblyContext context)
inlineprotectedinherited

Definition at line 67 of file elastic_cable_base.h.

References GRINS::ElasticCableBase::_disp_vars, and GRINS::SolidMechanicsVariables::u_var().

Referenced by GRINS::ElasticCableConstantGravity::element_time_derivative(), and GRINS::ElasticCableBase::init_context().

68  {
69  // For this Physics, we need to make sure that we grab only the 1D elements
70  return context.get_element_fe(_disp_vars.u_var(),1);
71  }
SolidMechanicsFEVariables _disp_vars
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.
void GRINS::ElasticCableBase::init_context ( AssemblyContext context)
virtualinherited

Initialize context for added physics variables.

Reimplemented from GRINS::Physics.

Definition at line 70 of file elastic_cable_base.C.

References GRINS::ElasticCableBase::get_fe().

71  {
72  this->get_fe(context)->get_JxW();
73  this->get_fe(context)->get_phi();
74  this->get_fe(context)->get_dphidxi();
75 
76  // Need for constructing metric tensors
77  this->get_fe(context)->get_dxyzdxi();
78  this->get_fe(context)->get_dxidx();
79  this->get_fe(context)->get_dxidy();
80  this->get_fe(context)->get_dxidz();
81 
82  return;
83  }
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)
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.
void GRINS::ElasticCableBase::init_variables ( libMesh::FEMSystem *  system)
virtualinherited

Initialize variables for this physics.

Implements GRINS::Physics.

Definition at line 51 of file elastic_cable_base.C.

References GRINS::ElasticCableBase::_disp_vars, and GRINS::SolidMechanicsFEVariables::init().

52  {
53  // is_2D = false, is_3D = true
54  _disp_vars.init(system,false,true);
55 
56  return;
57  }
SolidMechanicsFEVariables _disp_vars
void init(libMesh::FEMSystem *system, bool is_2D, bool is_3D)
Initialize System variables.
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
template<typename StressStrainLaw >
void GRINS::ElasticCable< StressStrainLaw >::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.

Definition at line 299 of file elastic_cable.C.

302  {
303  libmesh_not_implemented();
304  return;
305  }
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  }
void GRINS::Physics::read_input_options ( const GetPot &  input)
virtualinherited

Read options from GetPot input file. By default, nothing is read.

Reimplemented in GRINS::ScalarODE, GRINS::AveragedTurbineBase< Viscosity >, GRINS::AxisymmetricBoussinesqBuoyancy, GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::AxisymmetricHeatTransfer< Conductivity >, GRINS::AveragedFanBase< Viscosity >, GRINS::ParsedVelocitySourceBase< Viscosity >, GRINS::VelocityDragBase< Viscosity >, GRINS::VelocityPenaltyBase< Viscosity >, GRINS::IncompressibleNavierStokes< Viscosity >, GRINS::LowMachNavierStokes< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::HeatTransfer< Conductivity >, GRINS::ReactingLowMachNavierStokesBase< Mixture, Evaluator >, and GRINS::ReactingLowMachNavierStokes< Mixture, Evaluator >.

Definition at line 70 of file physics.C.

References GRINS::Physics::_enabled_subdomains, and GRINS::Physics::_physics_name.

Referenced by GRINS::HeatTransferBase< Conductivity >::HeatTransferBase(), GRINS::HeatTransferStabilizationBase< Conductivity >::HeatTransferStabilizationBase(), GRINS::IncompressibleNavierStokesAdjointStabilization< Viscosity >::IncompressibleNavierStokesAdjointStabilization(), GRINS::IncompressibleNavierStokesSPGSMStabilization< Viscosity >::IncompressibleNavierStokesSPGSMStabilization(), GRINS::Physics::Physics(), and GRINS::SpalartAllmarasSPGSMStabilization< Viscosity >::SpalartAllmarasSPGSMStabilization().

71  {
72  int num_ids = input.vector_variable_size( "Physics/"+this->_physics_name+"/enabled_subdomains" );
73 
74  for( int i = 0; i < num_ids; i++ )
75  {
76  libMesh::subdomain_id_type dumvar = input( "Physics/"+this->_physics_name+"/enabled_subdomains", -1, i );
77  _enabled_subdomains.insert( dumvar );
78  }
79 
80  return;
81  }
const PhysicsName _physics_name
Name of the physics object. Used for reading physics specific inputs.
Definition: physics.h:254
std::set< libMesh::subdomain_id_type > _enabled_subdomains
Subdomains on which the current Physics class is enabled.
Definition: physics.h:261
void GRINS::ParameterUser::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 in GRINS::AxisymmetricHeatTransfer< Conductivity >, GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >, GRINS::IncompressibleNavierStokesBase< Viscosity >, GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >, GRINS::HeatConduction< Conductivity >, GRINS::HeatTransferBase< Conductivity >, and GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >.

Definition at line 50 of file parameter_user.C.

Referenced by GRINS::BoussinesqBuoyancyAdjointStabilization< Viscosity >::register_parameter(), GRINS::HeatTransferBase< Conductivity >::register_parameter(), GRINS::HeatConduction< Conductivity >::register_parameter(), GRINS::BoussinesqBuoyancySPGSMStabilization< Viscosity >::register_parameter(), GRINS::IncompressibleNavierStokesBase< Viscosity >::register_parameter(), GRINS::LowMachNavierStokesBase< Viscosity, SpecificHeat, ThermalConductivity >::register_parameter(), and GRINS::AxisymmetricHeatTransfer< Conductivity >::register_parameter().

53  {
54  std::map<std::string, libMesh::Number*>::const_iterator it =
55  _my_parameters.find(param_name);
56 
57  if (it != _my_parameters.end())
58  {
59  std::cout << _my_name << " uses parameter " << param_name
60  << std::endl;
61  param_pointer.push_back(it->second);
62  }
63  }
std::map< std::string, libMesh::Number * > _my_parameters
template<typename StressStrainLaw >
void GRINS::ElasticCable< StressStrainLaw >::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtual

Register postprocessing variables for ElasticCable.

Reimplemented from GRINS::Physics.

Definition at line 79 of file elastic_cable.C.

References GRINS::elastic_cable, and GRINS::PostProcessedQuantities< NumericType >::register_quantity().

81  {
82  std::string section = "Physics/"+elastic_cable+"/output_vars";
83 
84  if( input.have_variable(section) )
85  {
86  unsigned int n_vars = input.vector_variable_size(section);
87 
88  for( unsigned int v = 0; v < n_vars; v++ )
89  {
90  std::string name = input(section,"DIE!",v);
91 
92  if( name == std::string("stress") )
93  {
94  // sigma_xx only, i.e., locally normal to the cutting plane
95  // sigma_yy=sigma_zz = 0 by assumption of this Physics
96  _stress_indices.resize(1);
97 
98  this->_stress_indices[0] = postprocessing.register_quantity("cable_stress");
99 
100  }
101  else if( name == std::string("strain") )
102  {
103  // eps_xx only, i.e., locally normal to the cutting plane
104  // eps_yy=eps_zz = 0 by assumption of this Physics
105  _strain_indices.resize(1);
106 
107  this->_strain_indices[0] = postprocessing.register_quantity("cable_strain");
108 
109  }
110  else if( name == std::string("force") )
111  {
112  // force_x only, i.e., locally normal to the cutting plane
113  // force_y=force_z=0 by assumption of this Physics
114  _force_indices.resize(1);
115 
116  this->_force_indices[0] = postprocessing.register_quantity("cable_force");
117 
118  }
119  else
120  {
121  std::cerr << "Error: Invalue output_vars value for "+elastic_cable << std::endl
122  << " Found " << name << std::endl
123  << " Acceptable values are: stress" << std::endl
124  << " strain" << std::endl
125  << " force " << std::endl;
126  libmesh_error();
127  }
128  }
129  }
130  return;
131  }
std::vector< unsigned int > _force_indices
Index from registering this quantity. Each component will have it's own index.
const PhysicsName elastic_cable
std::vector< unsigned int > _stress_indices
Index from registering this quantity. Each component will have it's own index.
Definition: elastic_cable.h:98
std::vector< unsigned int > _strain_indices
Index from registering this quantity. Each component will have it's own index.
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
void GRINS::ElasticCableBase::set_time_evolving_vars ( libMesh::FEMSystem *  system)
virtualinherited

Set which variables are time evolving.

Set those variables which evolve in time (as opposed to variables that behave like constraints). This is done separately from init_variables() because the MultiphysicsSystem must initialize its base class before time_evolving variables can be set. Default implementation is no time evolving variables.

Reimplemented from GRINS::Physics.

Definition at line 60 of file elastic_cable_base.C.

References GRINS::ElasticCableBase::_disp_vars, GRINS::SolidMechanicsVariables::u_var(), GRINS::SolidMechanicsVariables::v_var(), and GRINS::SolidMechanicsVariables::w_var().

61  {
62  // Tell the system to march temperature forward in time
63  system->time_evolving(_disp_vars.u_var());
64  system->time_evolving(_disp_vars.v_var());
65  system->time_evolving(_disp_vars.w_var());
66 
67  return;
68  }
SolidMechanicsFEVariables _disp_vars
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  }
template<typename StressStrainLaw >
void GRINS::ElasticCable< StressStrainLaw >::side_time_derivative ( bool  compute_jacobian,
AssemblyContext context,
CachedValues cache 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 281 of file elastic_cable.C.

284  {
285  std::vector<BoundaryID> ids = context.side_boundary_ids();
286 
287  for( std::vector<BoundaryID>::const_iterator it = ids.begin();
288  it != ids.end(); it++ )
289  {
290  libmesh_assert (*it != libMesh::BoundaryInfo::invalid_id);
291 
292  _bc_handler->apply_neumann_bcs( context, cache, compute_jacobian, *it );
293  }
294 
295  return;
296  }
GRINS::BCHandlingBase * _bc_handler
Definition: physics.h:256
virtual void apply_neumann_bcs(AssemblyContext &context, const GRINS::CachedValues &cache, const bool request_jacobian, const GRINS::BoundaryID bc_id) const

Member Data Documentation

template<typename StressStrainLaw >
libMesh::Real GRINS::ElasticCable< StressStrainLaw >::_A
private
GRINS::BCHandlingBase* GRINS::Physics::_bc_handler
protectedinherited
SolidMechanicsFEVariables GRINS::ElasticCableBase::_disp_vars
protectedinherited
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<typename StressStrainLaw >
std::vector<unsigned int> GRINS::ElasticCable< StressStrainLaw >::_force_indices
private

Index from registering this quantity. Each component will have it's own index.

Definition at line 104 of file elastic_cable.h.

GRINS::ICHandlingBase* GRINS::Physics::_ic_handler
protectedinherited
bool GRINS::Physics::_is_axisymmetric
protectedinherited
template<typename StressStrainLaw >
bool GRINS::ElasticCable< StressStrainLaw >::_is_compressible
private

Definition at line 95 of file elastic_cable.h.

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().

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<typename StressStrainLaw >
std::vector<unsigned int> GRINS::ElasticCable< StressStrainLaw >::_strain_indices
private

Index from registering this quantity. Each component will have it's own index.

Definition at line 101 of file elastic_cable.h.

template<typename StressStrainLaw >
std::vector<unsigned int> GRINS::ElasticCable< StressStrainLaw >::_stress_indices
private

Index from registering this quantity. Each component will have it's own index.

Definition at line 98 of file elastic_cable.h.

template<typename StressStrainLaw >
StressStrainLaw GRINS::ElasticCable< StressStrainLaw >::_stress_strain_law
private

Definition at line 91 of file elastic_cable.h.


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

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