GRINS-0.8.0
List of all members | Public Member Functions | Private Member Functions | Private Attributes
GRINS::ElasticMembrane< StressStrainLaw > Class Template Reference

#include <elastic_membrane.h>

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

Public Member Functions

 ElasticMembrane (const GRINS::PhysicsName &physics_name, const GetPot &input, bool is_compressible)
 
virtual ~ElasticMembrane ()
 
virtual void register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register postprocessing variables for ElasticMembrane. More...
 
virtual void element_time_derivative (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context)
 Constraint part(s) of physics for element interiors. More...
 
virtual void mass_residual (bool compute_jacobian, AssemblyContext &context)
 Mass matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. 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...
 
- Public Member Functions inherited from GRINS::ElasticMembraneBase< StressStrainLaw >
 ElasticMembraneBase (const GRINS::PhysicsName &physics_name, const GetPot &input, bool is_compressible)
 
virtual ~ElasticMembraneBase ()
 
virtual void init_variables (libMesh::FEMSystem *system)
 Initialize variables for this physics. More...
 
- Public Member Functions inherited from GRINS::ElasticMembraneAbstract
 ElasticMembraneAbstract (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~ElasticMembraneAbstract ()
 
virtual void init_context (AssemblyContext &context)
 Initialize context for added physics variables. More...
 
- Public Member Functions inherited from GRINS::SolidMechanicsAbstract
 SolidMechanicsAbstract (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~SolidMechanicsAbstract ()
 
virtual void set_time_evolving_vars (libMesh::FEMSystem *system)
 Set which variables are time evolving. More...
 
- 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 preassembly (MultiphysicsSystem &)
 Perform any necessary setup before element assembly begins. More...
 
virtual void reinit (MultiphysicsSystem &)
 Any reinitialization that needs to be done. More...
 
virtual void side_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void side_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool, AssemblyContext &)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void nonlocal_mass_residual (bool, AssemblyContext &)
 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 (AssemblyContext &)
 
virtual void compute_side_time_derivative_cache (AssemblyContext &)
 
virtual void compute_nonlocal_time_derivative_cache (AssemblyContext &)
 
virtual void compute_element_constraint_cache (AssemblyContext &)
 
virtual void compute_side_constraint_cache (AssemblyContext &)
 
virtual void compute_nonlocal_constraint_cache (AssemblyContext &)
 
virtual void compute_damping_residual_cache (AssemblyContext &)
 
virtual void compute_mass_residual_cache (AssemblyContext &)
 
virtual void compute_nonlocal_mass_residual_cache (AssemblyContext &)
 
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...
 
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...
 

Private Member Functions

 ElasticMembrane ()
 

Private Attributes

std::vector< unsigned int > _stress_indices
 Index from registering this quantity for postprocessing. Each component will have it's own index. More...
 
unsigned int _stress_zz_index
 Index from registering sigma_zz for postprocessing. Mainly for sanity checking. More...
 
std::vector< unsigned int > _strain_indices
 Index from registering this quantity for postprocessing. Each component will have it's own index. More...
 

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 Types inherited from GRINS::SolidMechanicsAbstract
typedef const libMesh::DenseSubVector< libMesh::Number > &(libMesh::DiffContext::* VarFuncType) (unsigned int) const
 
typedef void(libMesh::FEMContext::* InteriorFuncType) (unsigned int, unsigned int, libMesh::Real &) const
 
typedef libMesh::Real(libMesh::DiffContext::* VarDerivType) () const
 
- Protected Member Functions inherited from GRINS::ElasticMembraneBase< StressStrainLaw >
void mass_residual_impl (bool compute_jacobian, AssemblyContext &context, InteriorFuncType interior_solution, VarDerivType get_solution_deriv, libMesh::Real mu=1.0)
 Implementation of mass_residual. More...
 
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)
 
- Protected Member Functions inherited from GRINS::ElasticMembraneAbstract
const libMesh::FEGenericBase< libMesh::Real > * get_fe (const AssemblyContext &context)
 
- 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)
 
void check_var_subdomain_consistency (const FEVariablesBase &var) const
 Check that var is enabled on at least the subdomains this Physics is. More...
 
- Protected Attributes inherited from GRINS::ElasticMembraneBase< StressStrainLaw >
StressStrainLaw _stress_strain_law
 
bool _is_compressible
 
libMesh::Real _h0
 Membrane thickness. More...
 
libMesh::Real _rho
 Membrane density. More...
 
VariableIndex _lambda_sq_var
 Variable index for lambda_sq variable. More...
 
- Protected Attributes inherited from GRINS::SolidMechanicsAbstract
DisplacementVariable_disp_vars
 
- 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...
 
- 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<typename StressStrainLaw>
class GRINS::ElasticMembrane< StressStrainLaw >

Definition at line 34 of file elastic_membrane.h.

Constructor & Destructor Documentation

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

Definition at line 46 of file elastic_membrane.C.

References GRINS::Physics::_ic_handler.

48  : ElasticMembraneBase<StressStrainLaw>(physics_name,input,is_compressible)
49  {
50  this->_ic_handler = new GenericICHandler(physics_name, input);
51  }
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
template<typename StressStrainLaw >
virtual GRINS::ElasticMembrane< StressStrainLaw >::~ElasticMembrane ( )
inlinevirtual

Definition at line 41 of file elastic_membrane.h.

41 {};
template<typename StressStrainLaw >
GRINS::ElasticMembrane< StressStrainLaw >::ElasticMembrane ( )
private

Member Function Documentation

template<typename StressStrainLaw >
void GRINS::ElasticMembrane< 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 380 of file elastic_membrane.C.

References GRINS::Constants::two_pi.

384  {
385  bool is_stress = false;
386  if( !_stress_indices.empty() )
387  is_stress= ( _stress_indices[0] == quantity_index ||
388  _stress_indices[1] == quantity_index ||
389  _stress_indices[2] == quantity_index ||
390  _stress_indices[3] == quantity_index ||
391  _stress_indices[4] == quantity_index ||
392  _stress_indices[5] == quantity_index ||
393  _stress_indices[6] == quantity_index ||
394  _stress_zz_index == quantity_index );
395 
396  bool is_strain = false;
397  if( !_strain_indices.empty() )
398  is_strain = ( _strain_indices[0] == quantity_index ||
399  _strain_indices[1] == quantity_index ||
400  _strain_indices[2] == quantity_index );
401 
402  if( is_stress || is_strain )
403  {
404  const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
405 
406  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
407  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( this->_disp_vars.v() );
408  const libMesh::DenseSubVector<libMesh::Number>* w_coeffs = NULL;
409 
410  if( this->_disp_vars.dim() == 3 )
411  w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
412 
413  // Build new FE for the current point. We need this to build tensors at point.
414  libMesh::UniquePtr<libMesh::FEGenericBase<libMesh::Real> > fe_new =
415  this->build_new_fe( &context.get_elem(), this->get_fe(context),
416  point );
417 
418  const std::vector<std::vector<libMesh::Real> >& dphi_dxi =
419  fe_new->get_dphidxi();
420 
421  const std::vector<std::vector<libMesh::Real> >& dphi_deta =
422  fe_new->get_dphideta();
423 
424  // Need these to build up the covariant and contravariant metric tensors
425  const std::vector<libMesh::RealGradient>& dxdxi = fe_new->get_dxyzdxi();
426  const std::vector<libMesh::RealGradient>& dxdeta = fe_new->get_dxyzdeta();
427 
428  // Gradients are w.r.t. master element coordinates
429  libMesh::Gradient grad_u, grad_v, grad_w;
430  for( unsigned int d = 0; d < n_u_dofs; d++ )
431  {
432  libMesh::RealGradient u_gradphi( dphi_dxi[d][0], dphi_deta[d][0] );
433  grad_u += u_coeffs(d)*u_gradphi;
434  grad_v += v_coeffs(d)*u_gradphi;
435 
436  if( this->_disp_vars.dim() == 3 )
437  grad_w += (*w_coeffs)(d)*u_gradphi;
438  }
439 
440  libMesh::RealGradient grad_x( dxdxi[0](0), dxdeta[0](0) );
441  libMesh::RealGradient grad_y( dxdxi[0](1), dxdeta[0](1) );
442  libMesh::RealGradient grad_z( dxdxi[0](2), dxdeta[0](2) );
443 
444  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
445  libMesh::Real lambda_sq = 0;
446 
447  // We're only computing one point at a time, so qp = 0 always
448  this->compute_metric_tensors(0, *fe_new, context, grad_u, grad_v, grad_w,
449  a_cov, a_contra, A_cov, A_contra, lambda_sq );
450 
451  // We have everything we need for strain now, so check if we are computing strain
452  if( is_strain )
453  {
454  if( _strain_indices[0] == quantity_index )
455  {
456  value = 0.5*(A_cov(0,0) - a_cov(0,0));
457  }
458  else if( _strain_indices[1] == quantity_index )
459  {
460  value = 0.5*(A_cov(0,1) - a_cov(0,1));
461  }
462  else if( _strain_indices[2] == quantity_index )
463  {
464  value = 0.5*(A_cov(1,1) - a_cov(1,1));
465  }
466  else
467  {
468  //Wat?!
469  libmesh_error();
470  }
471  return;
472  }
473 
474  if( is_stress )
475  {
476  libMesh::Real det_a = a_cov(0,0)*a_cov(1,1) - a_cov(0,1)*a_cov(1,0);
477  libMesh::Real det_A = A_cov(0,0)*A_cov(1,1) - A_cov(0,1)*A_cov(1,0);
478 
479  libMesh::Real I3 = lambda_sq*det_A/det_a;
480 
482  this->_stress_strain_law.compute_stress(2,a_contra,a_cov,A_contra,A_cov,tau);
483 
484  if( _stress_indices[0] == quantity_index )
485  {
486  // Need to convert to Cauchy stress
487  value = tau(0,0)/std::sqrt(I3);
488  }
489  else if( _stress_indices[1] == quantity_index )
490  {
491  // Need to convert to Cauchy stress
492  value = tau(0,1)/std::sqrt(I3);
493  }
494  else if( _stress_indices[2] == quantity_index )
495  {
496  // Need to convert to Cauchy stress
497  value = tau(1,1)/std::sqrt(I3);
498  }
499  else if( _stress_indices[3] == quantity_index )
500  {
501  value = 0.5*(tau(0,0) + tau(1,1)) + std::sqrt(0.25*(tau(0,0)-tau(1,1))*(tau(0,0)-tau(1,1))
502  + tau(0,1)*tau(0,1) );
503  }
504  else if( _stress_indices[4] == quantity_index ||
505  _stress_indices[5] == quantity_index ||
506  _stress_indices[6] == quantity_index )
507  {
508  if(this->_is_compressible)
509  {
510  tau(2,2) = this->_stress_strain_law.compute_33_stress( a_contra, a_cov, A_contra, A_cov );
511  }
512 
513  libMesh::Real stress_I1 = tau(0,0) + tau(1,1) + tau(2,2);
514  libMesh::Real stress_I2 = 0.5*(stress_I1*stress_I1 - (tau(0,0)*tau(0,0) + tau(1,1)*tau(1,1)
515  + tau(2,2)*tau(2,2) + tau(0,1)*tau(0,1)
516  + tau(1,0)*tau(1,0)) );
517 
518  libMesh::Real stress_I3 = tau(2,2)*( tau(0,0)*tau(1,1) - tau(1,0)*tau(0,1) );
519 
520  /* Formulae for principal stresses from:
521  http://en.wikiversity.org/wiki/Principal_stresses */
522 
523  // I_2^2 - 3*I_2
524  libMesh::Real C1 = (stress_I1*stress_I1 - 3*stress_I2);
525 
526  // 2*I_1^3 - 9*I_1*_I2 + 27*I_3
527  libMesh::Real C2 = (2*stress_I1*stress_I1*stress_I1 - 9*stress_I1*stress_I2 + 27*stress_I3)/54;
528 
529  libMesh::Real theta = std::acos( C2/(2*std::sqrt(C1*C1*C1)) )/3.0;
530 
531  if( _stress_indices[4] == quantity_index )
532  {
533  value = (stress_I1 + 2.0*std::sqrt(C1)*std::cos(theta))/3.0;
534  }
535 
536  if( _stress_indices[5] == quantity_index )
537  {
538  value = (stress_I1 + 2.0*std::sqrt(C1)*std::cos(theta+Constants::two_pi/3.0))/3.0;
539  }
540 
541  if( _stress_indices[6] == quantity_index )
542  {
543  value = (stress_I1 + 2.0*std::sqrt(C1)*std::cos(theta+2.0*Constants::two_pi/3.0))/3.0;
544  }
545  }
546  else if( _stress_zz_index == quantity_index )
547  {
548  value = this->_stress_strain_law.compute_33_stress( a_contra, a_cov, A_contra, A_cov );
549  }
550  else
551  {
552  //Wat?!
553  libmesh_error();
554  }
555  } // is_stress
556 
557  }
558  }
unsigned int _stress_zz_index
Index from registering sigma_zz for postprocessing. Mainly for sanity checking.
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)
libMesh::UniquePtr< libMesh::FEGenericBase< libMesh::Real > > build_new_fe(const libMesh::Elem *elem, const libMesh::FEGenericBase< libMesh::Real > *fe, const libMesh::Point p)
Definition: physics.C:143
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)
std::vector< unsigned int > _stress_indices
Index from registering this quantity for postprocessing. Each component will have it's own index...
std::vector< unsigned int > _strain_indices
Index from registering this quantity for postprocessing. Each component will have it's own index...
const libMesh::Real two_pi
libMesh::Real compute_33_stress(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)
This is primarily a helper function for the plane stress cases.
unsigned int dim() const
Number of components.
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)
template<typename StressStrainLaw >
void GRINS::ElasticMembrane< StressStrainLaw >::element_constraint ( bool  ,
AssemblyContext  
)
virtual

Constraint part(s) of physics for element interiors.

Reimplemented from GRINS::Physics.

Definition at line 313 of file elastic_membrane.C.

314  {
315  // Only compute the constraint is tracking lambda_sq as an independent variable
316  if( this->_is_compressible )
317  {
318  unsigned int n_qpoints = context.get_element_qrule().n_points();
319 
320  const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
321 
322  const std::vector<libMesh::Real> &JxW = context.get_element_fe(this->_lambda_sq_var)->get_JxW();
323 
324  libMesh::DenseSubVector<libMesh::Number>& Fl = context.get_elem_residual(this->_lambda_sq_var);
325 
326  const std::vector<std::vector<libMesh::Real> >& phi =
327  context.get_element_fe(this->_lambda_sq_var)->get_phi();
328 
329  const unsigned int n_lambda_sq_dofs = context.get_dof_indices(this->_lambda_sq_var).size();
330 
331  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
332  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( this->_disp_vars.v() );
333  const libMesh::DenseSubVector<libMesh::Number>* w_coeffs = NULL;
334 
335  if( this->_disp_vars.dim() == 3 )
336  w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
337 
338  // All shape function gradients are w.r.t. master element coordinates
339  const std::vector<std::vector<libMesh::Real> >& dphi_dxi =
340  this->get_fe(context)->get_dphidxi();
341 
342  const std::vector<std::vector<libMesh::Real> >& dphi_deta =
343  this->get_fe(context)->get_dphideta();
344 
345  for (unsigned int qp=0; qp != n_qpoints; qp++)
346  {
347  libMesh::Real jac = JxW[qp];
348 
349  libMesh::Gradient grad_u, grad_v, grad_w;
350  for( unsigned int d = 0; d < n_u_dofs; d++ )
351  {
352  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp], dphi_deta[d][qp] );
353  grad_u += u_coeffs(d)*u_gradphi;
354  grad_v += v_coeffs(d)*u_gradphi;
355 
356  if( this->_disp_vars.dim() == 3 )
357  grad_w += (*w_coeffs)(d)*u_gradphi;
358  }
359 
360  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
361  libMesh::Real lambda_sq = 0;
362 
363  this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
364  grad_u, grad_v, grad_w,
365  a_cov, a_contra, A_cov, A_contra,
366  lambda_sq );
367 
368  libMesh::Real stress_33 = this->_stress_strain_law.compute_33_stress( a_contra, a_cov, A_contra, A_cov );
369 
370  for (unsigned int i=0; i != n_lambda_sq_dofs; i++)
371  Fl(i) += stress_33*phi[i][qp]*jac;
372 
373  if( compute_jacobian )
374  libmesh_not_implemented();
375  }
376  } // is_compressible
377  }
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)
VariableIndex _lambda_sq_var
Variable index for lambda_sq variable.
libMesh::Real compute_33_stress(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)
This is primarily a helper function for the plane stress cases.
unsigned int dim() const
Number of components.
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)
template<typename StressStrainLaw >
void GRINS::ElasticMembrane< StressStrainLaw >::element_time_derivative ( bool  compute_jacobian,
AssemblyContext context 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 117 of file elastic_membrane.C.

118  {
119  const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
120 
121  const std::vector<libMesh::Real> &JxW =
122  this->get_fe(context)->get_JxW();
123 
124  // Residuals that we're populating
125  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_disp_vars.u());
126  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_disp_vars.v());
127  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
128 
129  libMesh::DenseSubMatrix<libMesh::Number>& Kuu = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.u());
130  libMesh::DenseSubMatrix<libMesh::Number>& Kuv = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.v());
131 
132  libMesh::DenseSubMatrix<libMesh::Number>& Kvu = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.u());
133  libMesh::DenseSubMatrix<libMesh::Number>& Kvv = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.v());
134 
135  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
136  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
137 
138  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
139  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
140  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
141 
142  if( this->_disp_vars.dim() == 3 )
143  {
144  Fw = &context.get_elem_residual(this->_disp_vars.w());
145 
146  Kuw = &context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.w());
147  Kvw = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.w());
148  Kwu = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.u());
149  Kwv = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.v());
150  Kww = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.w());
151  }
152 
153  unsigned int n_qpoints = context.get_element_qrule().n_points();
154 
155  // All shape function gradients are w.r.t. master element coordinates
156  const std::vector<std::vector<libMesh::Real> >& dphi_dxi =
157  this->get_fe(context)->get_dphidxi();
158 
159  const std::vector<std::vector<libMesh::Real> >& dphi_deta =
160  this->get_fe(context)->get_dphideta();
161 
162  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
163  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( this->_disp_vars.v() );
164  const libMesh::DenseSubVector<libMesh::Number>* w_coeffs = NULL;
165 
166  if( this->_disp_vars.dim() == 3 )
167  w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
168 
169  // Need these to build up the covariant and contravariant metric tensors
170  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
171  const std::vector<libMesh::RealGradient>& dxdeta = this->get_fe(context)->get_dxyzdeta();
172 
173  for (unsigned int qp=0; qp != n_qpoints; qp++)
174  {
175  // Gradients are w.r.t. master element coordinates
176  libMesh::Gradient grad_u, grad_v, grad_w;
177 
178  for( unsigned int d = 0; d < n_u_dofs; d++ )
179  {
180  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp], dphi_deta[d][qp] );
181  grad_u += u_coeffs(d)*u_gradphi;
182  grad_v += v_coeffs(d)*u_gradphi;
183 
184  if( this->_disp_vars.dim() == 3 )
185  grad_w += (*w_coeffs)(d)*u_gradphi;
186  }
187 
188  libMesh::RealGradient grad_x( dxdxi[qp](0), dxdeta[qp](0) );
189  libMesh::RealGradient grad_y( dxdxi[qp](1), dxdeta[qp](1) );
190  libMesh::RealGradient grad_z( dxdxi[qp](2), dxdeta[qp](2) );
191 
192  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
193  libMesh::Real lambda_sq = 0;
194 
195  this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
196  grad_u, grad_v, grad_w,
197  a_cov, a_contra, A_cov, A_contra,
198  lambda_sq );
199 
200  const unsigned int manifold_dim = 2; // The manifold dimension is always 2 for this physics
201 
202  // Compute stress and elasticity tensors
204  ElasticityTensor C;
205  this->_stress_strain_law.compute_stress_and_elasticity(manifold_dim,a_contra,a_cov,A_contra,A_cov,tau,C);
206 
207  libMesh::Real jac = JxW[qp];
208 
209  for (unsigned int i=0; i != n_u_dofs; i++)
210  {
211  libMesh::RealGradient u_gradphi( dphi_dxi[i][qp], dphi_deta[i][qp] );
212 
213  for( unsigned int alpha = 0; alpha < manifold_dim; alpha++ )
214  {
215  for( unsigned int beta = 0; beta < manifold_dim; beta++ )
216  {
217  libMesh::Real factor = 0.5*tau(alpha,beta)*this->_h0*jac;
218 
219  Fu(i) += factor*( (grad_x(beta) + grad_u(beta))*u_gradphi(alpha) +
220  (grad_x(alpha) + grad_u(alpha))*u_gradphi(beta) );
221 
222  Fv(i) += factor*( (grad_y(beta) + grad_v(beta))*u_gradphi(alpha) +
223  (grad_y(alpha) + grad_v(alpha))*u_gradphi(beta) );
224 
225  if( this->_disp_vars.dim() == 3 )
226  (*Fw)(i) += factor*( (grad_z(beta) + grad_w(beta))*u_gradphi(alpha) +
227  (grad_z(alpha) + grad_w(alpha))*u_gradphi(beta) );
228  }
229  }
230  }
231 
232  if( compute_jacobian )
233  {
234  for (unsigned int i=0; i != n_u_dofs; i++)
235  {
236  libMesh::RealGradient u_gradphi_i( dphi_dxi[i][qp], dphi_deta[i][qp] );
237 
238  for (unsigned int j=0; j != n_u_dofs; j++)
239  {
240  libMesh::RealGradient u_gradphi_j( dphi_dxi[j][qp], dphi_deta[j][qp] );
241 
242  for( unsigned int alpha = 0; alpha < manifold_dim; alpha++ )
243  {
244  for( unsigned int beta = 0; beta < manifold_dim; beta++ )
245  {
246  const libMesh::Real diag_term = 0.5*this->_h0*jac*tau(alpha,beta)*context.get_elem_solution_derivative()*
247  ( u_gradphi_j(beta)*u_gradphi_i(alpha) +
248  u_gradphi_j(alpha)*u_gradphi_i(beta) );
249  Kuu(i,j) += diag_term;
250 
251  Kvv(i,j) += diag_term;
252 
253  if( this->_disp_vars.dim() == 3 )
254  (*Kww)(i,j) += diag_term;
255 
256  for( unsigned int lambda = 0; lambda < manifold_dim; lambda++ )
257  {
258  for( unsigned int mu = 0; mu < manifold_dim; mu++ )
259  {
260  const libMesh::Real dgamma_du = 0.5*( u_gradphi_j(lambda)*(grad_x(mu)+grad_u(mu)) +
261  (grad_x(lambda)+grad_u(lambda))*u_gradphi_j(mu) );
262 
263  const libMesh::Real dgamma_dv = 0.5*( u_gradphi_j(lambda)*(grad_y(mu)+grad_v(mu)) +
264  (grad_y(lambda)+grad_v(lambda))*u_gradphi_j(mu) );
265 
266  const libMesh::Real C1 = 0.5*this->_h0*jac*C(alpha,beta,lambda,mu)*context.get_elem_solution_derivative();
267 
268  const libMesh::Real x_term = C1*( (grad_x(beta)+grad_u(beta))*u_gradphi_i(alpha) +
269  (grad_x(alpha)+grad_u(alpha))*u_gradphi_i(beta) );
270 
271  const libMesh::Real y_term = C1*( (grad_y(beta)+grad_v(beta))*u_gradphi_i(alpha) +
272  (grad_y(alpha)+grad_v(alpha))*u_gradphi_i(beta) );
273 
274  Kuu(i,j) += x_term*dgamma_du;
275 
276  Kuv(i,j) += x_term*dgamma_dv;
277 
278  Kvu(i,j) += y_term*dgamma_du;
279 
280  Kvv(i,j) += y_term*dgamma_dv;
281 
282  if( this->_disp_vars.dim() == 3 )
283  {
284  const libMesh::Real dgamma_dw = 0.5*( u_gradphi_j(lambda)*(grad_z(mu)+grad_w(mu)) +
285  (grad_z(lambda)+grad_w(lambda))*u_gradphi_j(mu) );
286 
287  const libMesh::Real z_term = C1*( (grad_z(beta)+grad_w(beta))*u_gradphi_i(alpha) +
288  (grad_z(alpha)+grad_w(alpha))*u_gradphi_i(beta) );
289 
290  (*Kuw)(i,j) += x_term*dgamma_dw;
291 
292  (*Kvw)(i,j) += y_term*dgamma_dw;
293 
294  (*Kwu)(i,j) += z_term*dgamma_du;
295 
296  (*Kwv)(i,j) += z_term*dgamma_dv;
297 
298  (*Kww)(i,j) += z_term*dgamma_dw;
299  }
300  }
301  }
302  }
303  }
304  }
305  }
306  }
307 
308  }
309  }
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_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)
unsigned int dim() const
Number of components.
libMesh::Real _h0
Membrane thickness.
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)
template<typename StressStrainLaw >
virtual void GRINS::ElasticMembrane< StressStrainLaw >::mass_residual ( bool  ,
AssemblyContext  
)
inlinevirtual

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

Reimplemented from GRINS::Physics.

Definition at line 54 of file elastic_membrane.h.

References GRINS::ElasticMembraneBase< StressStrainLaw >::mass_residual_impl().

56  { this->mass_residual_impl(compute_jacobian,
57  context,
58  &libMesh::FEMContext::interior_accel,
59  &libMesh::DiffContext::get_elem_solution_accel_derivative); }
void mass_residual_impl(bool compute_jacobian, AssemblyContext &context, InteriorFuncType interior_solution, VarDerivType get_solution_deriv, libMesh::Real mu=1.0)
Implementation of mass_residual.
template<typename StressStrainLaw >
void GRINS::ElasticMembrane< StressStrainLaw >::register_postprocessing_vars ( const GetPot &  input,
PostProcessedQuantities< libMesh::Real > &  postprocessing 
)
virtual

Register postprocessing variables for ElasticMembrane.

Reimplemented from GRINS::Physics.

Definition at line 54 of file elastic_membrane.C.

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

56  {
57  std::string section = "Physics/"+PhysicsNaming::elastic_membrane()+"/output_vars";
58 
59  if( input.have_variable(section) )
60  {
61  unsigned int n_vars = input.vector_variable_size(section);
62 
63  for( unsigned int v = 0; v < n_vars; v++ )
64  {
65  std::string name = input(section,"DIE!",v);
66 
67  if( name == std::string("stress") )
68  {
69  // sigma_xx, sigma_xy, sigma_yy, sigma_yx = sigma_xy
70  // sigma_zz = 0 by assumption of this Physics
71  _stress_indices.resize(7);
72 
73  this->_stress_indices[0] = postprocessing.register_quantity("stress_xx");
74 
75  this->_stress_indices[1] = postprocessing.register_quantity("stress_xy");
76 
77  this->_stress_indices[2] = postprocessing.register_quantity("stress_yy");
78 
79  this->_stress_indices[3] = postprocessing.register_quantity("sigma_max");
80 
81  this->_stress_indices[4] = postprocessing.register_quantity("sigma_1");
82 
83  this->_stress_indices[5] = postprocessing.register_quantity("sigma_2");
84 
85  this->_stress_indices[6] = postprocessing.register_quantity("sigma_3");
86  }
87  else if( name == std::string( "stress_zz" ) )
88  {
89  // This is mostly for sanity checking the plane stress condition
90  this->_stress_zz_index = postprocessing.register_quantity("stress_zz");
91  }
92  else if( name == std::string("strain") )
93  {
94  // eps_xx, eps_xy, eps_yy, eps_yx = eps_xy
95  _strain_indices.resize(3);
96 
97  this->_strain_indices[0] = postprocessing.register_quantity("strain_xx");
98 
99  this->_strain_indices[1] = postprocessing.register_quantity("strain_xy");
100 
101  this->_strain_indices[2] = postprocessing.register_quantity("strain_yy");
102  }
103  else
104  {
105  std::cerr << "Error: Invalue output_vars value for "+PhysicsNaming::elastic_membrane() << std::endl
106  << " Found " << name << std::endl
107  << " Acceptable values are: stress" << std::endl
108  << " strain" << std::endl;
109  libmesh_error();
110  }
111  }
112  }
113  }
unsigned int _stress_zz_index
Index from registering sigma_zz for postprocessing. Mainly for sanity checking.
std::vector< unsigned int > _stress_indices
Index from registering this quantity for postprocessing. Each component will have it's own index...
static PhysicsName elastic_membrane()
std::vector< unsigned int > _strain_indices
Index from registering this quantity for postprocessing. Each component will have it's own index...

Member Data Documentation

template<typename StressStrainLaw >
std::vector<unsigned int> GRINS::ElasticMembrane< StressStrainLaw >::_strain_indices
private

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

Definition at line 78 of file elastic_membrane.h.

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

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

Definition at line 72 of file elastic_membrane.h.

template<typename StressStrainLaw >
unsigned int GRINS::ElasticMembrane< StressStrainLaw >::_stress_zz_index
private

Index from registering sigma_zz for postprocessing. Mainly for sanity checking.

Definition at line 75 of file elastic_membrane.h.


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

Generated on Tue Dec 19 2017 12:47:30 for GRINS-0.8.0 by  doxygen 1.8.9.1