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

#include <elastic_membrane_base.h>

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

Public Member Functions

 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 register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. 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 nonlocal_time_derivative (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Time dependent part(s) of physics for scalar variables. More...
 
virtual void element_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for element interiors. More...
 
virtual void side_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for boundaries of elements on the domain boundary. More...
 
virtual void nonlocal_constraint (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Constraint part(s) of physics for scalar variables. More...
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
 Damping matrix part(s) for element interiors. All boundary terms lie within the time_derivative part. More...
 
virtual void 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_ics (libMesh::FEMSystem *system, libMesh::CompositeFunction< libMesh::Number > &all_ics)
 
virtual void compute_element_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_time_derivative_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_element_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_side_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_constraint_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_damping_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_nonlocal_mass_residual_cache (const AssemblyContext &context, CachedValues &cache)
 
virtual void compute_postprocessed_quantity (unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
 
ICHandlingBaseget_ic_handler ()
 
- Public Member Functions inherited from GRINS::ParameterUser
 ParameterUser (const std::string &user_name)
 
virtual ~ParameterUser ()
 
virtual void set_parameter (libMesh::Number &param_variable, const GetPot &input, const std::string &param_name, libMesh::Number param_default)
 Each subclass can simultaneously read a parameter value from. More...
 
virtual void set_parameter (libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void set_parameter (libMesh::ParsedFEMFunction< libMesh::Number > &func, const GetPot &input, const std::string &func_param_name, const std::string &param_default)
 Each subclass can simultaneously read a parsed function from. More...
 
virtual void move_parameter (const libMesh::Number &old_parameter, libMesh::Number &new_parameter)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &old_func, libMesh::ParsedFunction< libMesh::Number, libMesh::Gradient > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
virtual void move_parameter (const libMesh::ParsedFEMFunction< libMesh::Number > &old_func, libMesh::ParsedFEMFunction< libMesh::Number > &new_func)
 When cloning an object, we need to update parameter pointers. More...
 
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...
 

Protected Member Functions

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)
 

Protected Attributes

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
DisplacementFEVariables _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...
 

Private Member Functions

 ElasticMembraneBase ()
 

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
 
- 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::ElasticMembraneBase< StressStrainLaw >

Definition at line 34 of file elastic_membrane_base.h.

Constructor & Destructor Documentation

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

Definition at line 40 of file elastic_membrane_base.C.

References GRINS::ElasticMembraneBase< StressStrainLaw >::_h0, GRINS::ElasticMembraneBase< StressStrainLaw >::_rho, GRINS::PhysicsNaming::elastic_membrane(), GRINS::MaterialsParsing::read_density(), and GRINS::MaterialsParsing::read_property().

43  : ElasticMembraneAbstract(physics_name,input),
45  _is_compressible(is_compressible),
46  _h0(0.0)
47  {
49  "Physics/"+physics_name+"/h0",
50  "MembraneThickness",
52  (*this),
53  _h0 );
54 
55  MaterialsParsing::read_density( physics_name,
56  input,
57  (*this),
58  _rho );
59  }
libMesh::Real _rho
Membrane density.
static void read_density(const std::string &core_physics_name, const GetPot &input, ParameterUser &params, libMesh::Real &rho)
Helper function to reading density from input.
static PhysicsName elastic_membrane()
static void read_property(const GetPot &input, const std::string &old_option, const std::string &property, const std::string &core_physics, ParameterUser &param_user, libMesh::Real &value)
Helper function for parsing/maintaing backward compatibility.
static std::string material_name(const GetPot &input, const std::string &physics)
Get the name of the material in the Physics/physics section.
libMesh::Real _h0
Membrane thickness.
template<typename StressStrainLaw >
virtual GRINS::ElasticMembraneBase< StressStrainLaw >::~ElasticMembraneBase ( )
inlinevirtual

Definition at line 42 of file elastic_membrane_base.h.

42 {};
template<typename StressStrainLaw >
GRINS::ElasticMembraneBase< StressStrainLaw >::ElasticMembraneBase ( )
private

Member Function Documentation

template<typename StressStrainLaw >
void GRINS::ElasticMembraneBase< 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 
)
protected

Definition at line 133 of file elastic_membrane_base.C.

144  {
145  const std::vector<libMesh::RealGradient>& dxdxi = elem.get_dxyzdxi();
146  const std::vector<libMesh::RealGradient>& dxdeta = elem.get_dxyzdeta();
147 
148  const std::vector<libMesh::Real>& dxidx = elem.get_dxidx();
149  const std::vector<libMesh::Real>& dxidy = elem.get_dxidy();
150  const std::vector<libMesh::Real>& dxidz = elem.get_dxidz();
151 
152  const std::vector<libMesh::Real>& detadx = elem.get_detadx();
153  const std::vector<libMesh::Real>& detady = elem.get_detady();
154  const std::vector<libMesh::Real>& detadz = elem.get_detadz();
155 
156  libMesh::RealGradient dxi( dxidx[qp], dxidy[qp], dxidz[qp] );
157  libMesh::RealGradient deta( detadx[qp], detady[qp], detadz[qp] );
158 
159  libMesh::RealGradient dudxi( grad_u(0), grad_v(0), grad_w(0) );
160  libMesh::RealGradient dudeta( grad_u(1), grad_v(1), grad_w(1) );
161 
162  // Covariant metric tensor of reference configuration
163  a_cov.zero();
164  a_cov(0,0) = dxdxi[qp]*dxdxi[qp];
165  a_cov(0,1) = dxdxi[qp]*dxdeta[qp];
166  a_cov(1,0) = dxdeta[qp]*dxdxi[qp];
167  a_cov(1,1) = dxdeta[qp]*dxdeta[qp];
168 
169  libMesh::Real det_a = a_cov(0,0)*a_cov(1,1) - a_cov(0,1)*a_cov(1,0);
170 
171  // Covariant metric tensor of current configuration
172  A_cov.zero();
173  A_cov(0,0) = (dxdxi[qp] + dudxi)*(dxdxi[qp] + dudxi);
174  A_cov(0,1) = (dxdxi[qp] + dudxi)*(dxdeta[qp] + dudeta);
175  A_cov(1,0) = (dxdeta[qp] + dudeta)*(dxdxi[qp] + dudxi);
176  A_cov(1,1) = (dxdeta[qp] + dudeta)*(dxdeta[qp] + dudeta);
177 
178  // Contravariant metric tensor of reference configuration
179  a_contra.zero();
180  a_contra(0,0) = dxi*dxi;
181  a_contra(0,1) = dxi*deta;
182  a_contra(1,0) = deta*dxi;
183  a_contra(1,1) = deta*deta;
184 
185  // Contravariant metric tensor in current configuration is A_cov^{-1}
186  libMesh::Real det_A = A_cov(0,0)*A_cov(1,1) - A_cov(0,1)*A_cov(1,0);
187 
188  A_contra.zero();
189  A_contra(0,0) = A_cov(1,1)/det_A;
190  A_contra(0,1) = -A_cov(0,1)/det_A;
191  A_contra(1,0) = -A_cov(1,0)/det_A;
192  A_contra(1,1) = A_cov(0,0)/det_A;
193 
194  a_cov(2,2) = 1.0;
195  a_contra(2,2) = 1.0;
196 
197 
198  // If the material is compressible, then lambda_sq is an independent variable
199  if( _is_compressible )
200  {
201  lambda_sq = context.interior_value(this->_lambda_sq_var, qp);
202  }
203  else
204  {
205  // If the material is incompressible, lambda^2 is known
206  lambda_sq = det_a/det_A;
207  }
208 
209  A_cov(2,2) = lambda_sq;
210  A_contra(2,2) = 1.0/lambda_sq;
211  }
VariableIndex _lambda_sq_var
Variable index for lambda_sq variable.
template<typename StressStrainLaw >
void GRINS::ElasticMembraneBase< StressStrainLaw >::init_variables ( libMesh::FEMSystem *  system)
virtual

Initialize variables for this physics.

Todo:
Might want to make Order/FEType inputable

Reimplemented from GRINS::SolidMechanicsAbstract.

Definition at line 62 of file elastic_membrane_base.C.

References GRINS::SolidMechanicsAbstract::init_variables().

63  {
64  // First call base class
66 
67  // Now build lambda_sq variable if we need it
69  {
71  _lambda_sq_var = system->add_variable( "lambda_sq", GRINSEnums::FIRST, GRINSEnums::LAGRANGE);
72  }
73  }
virtual void init_variables(libMesh::FEMSystem *system)
Initialize variables for this physics.
VariableIndex _lambda_sq_var
Variable index for lambda_sq variable.
template<typename StressStrainLaw >
void GRINS::ElasticMembraneBase< StressStrainLaw >::mass_residual_impl ( bool  compute_jacobian,
AssemblyContext context,
InteriorFuncType  interior_solution,
VarDerivType  get_solution_deriv,
libMesh::Real  mu = 1.0 
)
protected

Implementation of mass_residual.

The mu argument is needed to support Rayleigh Damping.

Definition at line 76 of file elastic_membrane_base.C.

Referenced by GRINS::ElasticMembrane< StressStrainLaw >::mass_residual().

81  {
82  const unsigned int n_u_dofs = context.get_dof_indices(_disp_vars.u()).size();
83 
84  const std::vector<libMesh::Real> &JxW =
85  this->get_fe(context)->get_JxW();
86 
87  const std::vector<std::vector<libMesh::Real> >& u_phi =
88  this->get_fe(context)->get_phi();
89 
90  // Residuals that we're populating
91  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_disp_vars.u());
92  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_disp_vars.v());
93  libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(_disp_vars.w());
94 
95  libMesh::DenseSubMatrix<libMesh::Number>& Kuu = context.get_elem_jacobian(_disp_vars.u(),_disp_vars.u());
96  libMesh::DenseSubMatrix<libMesh::Number>& Kvv = context.get_elem_jacobian(_disp_vars.v(),_disp_vars.v());
97  libMesh::DenseSubMatrix<libMesh::Number>& Kww = context.get_elem_jacobian(_disp_vars.w(),_disp_vars.w());
98 
99  unsigned int n_qpoints = context.get_element_qrule().n_points();
100 
101  for (unsigned int qp=0; qp != n_qpoints; qp++)
102  {
103  libMesh::Real jac = JxW[qp];
104 
105  libMesh::Real u_ddot, v_ddot, w_ddot;
106  (context.*interior_solution)( _disp_vars.u(), qp, u_ddot );
107  (context.*interior_solution)( _disp_vars.v(), qp, v_ddot );
108  (context.*interior_solution)( _disp_vars.w(), qp, w_ddot );
109 
110  for (unsigned int i=0; i != n_u_dofs; i++)
111  {
112  Fu(i) += mu*this->_rho*_h0*u_ddot*u_phi[i][qp]*jac;
113  Fv(i) += mu*this->_rho*_h0*v_ddot*u_phi[i][qp]*jac;
114  Fw(i) += mu*this->_rho*_h0*w_ddot*u_phi[i][qp]*jac;
115 
116  if( compute_jacobian )
117  {
118  for (unsigned int j=0; j != n_u_dofs; j++)
119  {
120  libMesh::Real jac_term = this->_rho*_h0*u_phi[i][qp]*u_phi[j][qp]*jac;
121  jac_term *= mu*(context.*get_solution_deriv)();
122 
123  Kuu(i,j) += jac_term;
124  Kvv(i,j) += jac_term;
125  Kww(i,j) += jac_term;
126  }
127  }
128  }
129  }
130  }
libMesh::Real _rho
Membrane density.
libMesh::Real _h0
Membrane thickness.
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)

Member Data Documentation

template<typename StressStrainLaw >
libMesh::Real GRINS::ElasticMembraneBase< StressStrainLaw >::_h0
protected

Membrane thickness.

Definition at line 73 of file elastic_membrane_base.h.

Referenced by GRINS::ElasticMembraneBase< StressStrainLaw >::ElasticMembraneBase().

template<typename StressStrainLaw >
bool GRINS::ElasticMembraneBase< StressStrainLaw >::_is_compressible
protected

Definition at line 70 of file elastic_membrane_base.h.

template<typename StressStrainLaw >
VariableIndex GRINS::ElasticMembraneBase< StressStrainLaw >::_lambda_sq_var
protected

Variable index for lambda_sq variable.

Definition at line 79 of file elastic_membrane_base.h.

template<typename StressStrainLaw >
libMesh::Real GRINS::ElasticMembraneBase< StressStrainLaw >::_rho
protected

Membrane density.

Definition at line 76 of file elastic_membrane_base.h.

Referenced by GRINS::ElasticMembraneBase< StressStrainLaw >::ElasticMembraneBase().

template<typename StressStrainLaw >
StressStrainLaw GRINS::ElasticMembraneBase< StressStrainLaw >::_stress_strain_law
protected

Definition at line 68 of file elastic_membrane_base.h.


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

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