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

#include <elastic_membrane_rayleigh_damping.h>

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

Public Member Functions

 ElasticMembraneRayleighDamping (const PhysicsName &physics_name, const GetPot &input, bool is_compressible)
 
virtual ~ElasticMembraneRayleighDamping ()
 
virtual void damping_residual (bool compute_jacobian, AssemblyContext &context)
 Time dependent part(s) of physics for element interiors. 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 register_postprocessing_vars (const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
 Register name of postprocessed quantity with PostProcessedQuantities. 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 element_time_derivative (bool, AssemblyContext &)
 Time dependent part(s) of physics for element interiors. 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 element_constraint (bool, AssemblyContext &)
 Constraint part(s) of physics for element interiors. 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 mass_residual (bool, AssemblyContext &)
 Mass 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 &)
 
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 Attributes

libMesh::Real _lambda_factor
 
libMesh::Real _mu_factor
 
- 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...
 

Private Member Functions

 ElasticMembraneRayleighDamping ()
 

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...
 
- 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::ElasticMembraneRayleighDamping< StressStrainLaw >

Definition at line 33 of file elastic_membrane_rayleigh_damping.h.

Constructor & Destructor Documentation

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

Definition at line 40 of file elastic_membrane_rayleigh_damping.C.

References GRINS::PhysicsNaming::elastic_membrane(), GRINS::PhysicsNaming::elastic_membrane_rayleigh_damping(), and GRINS::Physics::parse_enabled_subdomains().

43  : ElasticMembraneBase<StressStrainLaw>(physics_name,input,is_compressible),
44  _lambda_factor(input("Physics/"+PhysicsNaming::elastic_membrane_rayleigh_damping()+"/lambda_factor",0.0)),
45  _mu_factor(input("Physics/"+PhysicsNaming::elastic_membrane_rayleigh_damping()+"/mu_factor",0.0))
46  {
47  if( !input.have_variable("Physics/"+PhysicsNaming::elastic_membrane_rayleigh_damping()+"/lambda_factor") )
48  libmesh_error_msg("ERROR: Couldn't find Physics/"+PhysicsNaming::elastic_membrane_rayleigh_damping()+"/lambda_factor in input!");
49 
50  if( !input.have_variable("Physics/"+PhysicsNaming::elastic_membrane_rayleigh_damping()+"/mu_factor") )
51  libmesh_error_msg("ERROR: Couldn't find Physics/"+PhysicsNaming::elastic_membrane_rayleigh_damping()+"/mu_factor in input!");
52 
53  // If the user specified enabled subdomains in this Physics section,
54  // that's an error; we're slave to ElasticMembrane.
55  if( input.have_variable("Physics/"+PhysicsNaming::elastic_membrane_rayleigh_damping()+"/enabled_subdomains" ) )
56  libmesh_error_msg("ERROR: Cannot specify subdomains for "
58  +"! Must specify subdomains through "
60 
62  }
static PhysicsName elastic_membrane()
void parse_enabled_subdomains(const GetPot &input, const std::string &physics_name)
Definition: physics.C:65
static PhysicsName elastic_membrane_rayleigh_damping()
template<typename StressStrainLaw >
virtual GRINS::ElasticMembraneRayleighDamping< StressStrainLaw >::~ElasticMembraneRayleighDamping ( )
inlinevirtual

Definition at line 41 of file elastic_membrane_rayleigh_damping.h.

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

Member Function Documentation

template<typename StressStrainLaw >
void GRINS::ElasticMembraneRayleighDamping< StressStrainLaw >::damping_residual ( bool  compute_jacobian,
AssemblyContext context 
)
virtual

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

Reimplemented from GRINS::Physics.

Definition at line 66 of file elastic_membrane_rayleigh_damping.C.

67  {
68  // First do the mass part
69  this->mass_residual_impl(compute_jacobian,
70  context,
71  &libMesh::FEMContext::interior_rate,
72  &libMesh::DiffContext::get_elem_solution_rate_derivative,
73  _mu_factor);
74 
75  // Now do the stiffness part
76  const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
77 
78  const std::vector<libMesh::Real> &JxW =
79  this->get_fe(context)->get_JxW();
80 
81  // Residuals that we're populating
82  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_disp_vars.u());
83  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_disp_vars.v());
84  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
85 
86  libMesh::DenseSubMatrix<libMesh::Number>& Kuu = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.u());
87  libMesh::DenseSubMatrix<libMesh::Number>& Kuv = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.v());
88 
89  libMesh::DenseSubMatrix<libMesh::Number>& Kvu = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.u());
90  libMesh::DenseSubMatrix<libMesh::Number>& Kvv = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.v());
91 
92  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
93  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
94 
95  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
96  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
97  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
98 
99  if( this->_disp_vars.dim() == 3 )
100  {
101  Fw = &context.get_elem_residual(this->_disp_vars.w());
102 
103  Kuw = &context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.w());
104  Kvw = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.w());
105  Kwu = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.u());
106  Kwv = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.v());
107  Kww = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.w());
108  }
109 
110  unsigned int n_qpoints = context.get_element_qrule().n_points();
111 
112  // All shape function gradients are w.r.t. master element coordinates
113  const std::vector<std::vector<libMesh::Real> >& dphi_dxi =
114  this->get_fe(context)->get_dphidxi();
115 
116  const std::vector<std::vector<libMesh::Real> >& dphi_deta =
117  this->get_fe(context)->get_dphideta();
118 
119  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
120  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( this->_disp_vars.v() );
121  const libMesh::DenseSubVector<libMesh::Number>* w_coeffs = NULL;
122 
123  const libMesh::DenseSubVector<libMesh::Number>& dudt_coeffs = context.get_elem_solution_rate( this->_disp_vars.u() );
124  const libMesh::DenseSubVector<libMesh::Number>& dvdt_coeffs = context.get_elem_solution_rate( this->_disp_vars.v() );
125  const libMesh::DenseSubVector<libMesh::Number>* dwdt_coeffs = NULL;
126 
127  if( this->_disp_vars.dim() == 3 )
128  {
129  w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
130  dwdt_coeffs = &context.get_elem_solution_rate( this->_disp_vars.w() );
131  }
132 
133  // Need these to build up the covariant and contravariant metric tensors
134  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
135  const std::vector<libMesh::RealGradient>& dxdeta = this->get_fe(context)->get_dxyzdeta();
136 
137  for (unsigned int qp=0; qp != n_qpoints; qp++)
138  {
139  // Gradients are w.r.t. master element coordinates
140  libMesh::Gradient grad_u, grad_v, grad_w;
141  libMesh::Gradient dgradu_dt, dgradv_dt, dgradw_dt;
142 
143  for( unsigned int d = 0; d < n_u_dofs; d++ )
144  {
145  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp], dphi_deta[d][qp] );
146  grad_u += u_coeffs(d)*u_gradphi;
147  grad_v += v_coeffs(d)*u_gradphi;
148 
149  dgradu_dt += dudt_coeffs(d)*u_gradphi;
150  dgradv_dt += dvdt_coeffs(d)*u_gradphi;
151 
152  if( this->_disp_vars.dim() == 3 )
153  {
154  grad_w += (*w_coeffs)(d)*u_gradphi;
155  dgradw_dt += (*dwdt_coeffs)(d)*u_gradphi;
156  }
157  }
158 
159  libMesh::RealGradient grad_x( dxdxi[qp](0), dxdeta[qp](0) );
160  libMesh::RealGradient grad_y( dxdxi[qp](1), dxdeta[qp](1) );
161  libMesh::RealGradient grad_z( dxdxi[qp](2), dxdeta[qp](2) );
162 
163 
164  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
165  libMesh::Real lambda_sq = 0;
166 
167  this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
168  grad_u, grad_v, grad_w,
169  a_cov, a_contra, A_cov, A_contra,
170  lambda_sq );
171 
172  const unsigned int manifold_dim = 2; // The manifold dimension is always 2 for this physics
173 
174  // Compute stress and elasticity tensors
176  ElasticityTensor C;
177  this->_stress_strain_law.compute_stress_and_elasticity(manifold_dim,a_contra,a_cov,A_contra,A_cov,tau,C);
178 
179  libMesh::Real jac = JxW[qp];
180 
181  for (unsigned int i=0; i != n_u_dofs; i++)
182  {
183  libMesh::RealGradient u_gradphi( dphi_dxi[i][qp], dphi_deta[i][qp] );
184 
185  for( unsigned int alpha = 0; alpha < manifold_dim; alpha++ )
186  {
187  for( unsigned int beta = 0; beta < manifold_dim; beta++ )
188  {
189  const libMesh::Real common_factor = _lambda_factor*0.5*this->_h0*jac;
190  const libMesh::Real factor = common_factor*tau(alpha,beta);
191 
192  libMesh::Real u_diagterm = dgradu_dt(beta)*u_gradphi(alpha) + dgradu_dt(alpha)*u_gradphi(beta);
193  libMesh::Real v_diagterm = dgradv_dt(beta)*u_gradphi(alpha) + dgradv_dt(alpha)*u_gradphi(beta);
194  libMesh::Real w_diagterm = dgradw_dt(beta)*u_gradphi(alpha) + dgradw_dt(alpha)*u_gradphi(beta);
195 
196  Fu(i) += factor*u_diagterm;
197  Fv(i) += factor*v_diagterm;
198 
199  if( this->_disp_vars.dim() == 3 )
200  (*Fw)(i) += factor*w_diagterm;
201 
202  for( unsigned int lambda = 0; lambda < manifold_dim; lambda++ )
203  {
204  for( unsigned int mu = 0; mu < manifold_dim; mu++ )
205  {
206  const libMesh::Real C1 = common_factor*C(alpha,beta,lambda,mu);
207 
208  const libMesh::Real gamma_u = 0.5*( dgradu_dt(lambda)*(grad_x(mu)+grad_u(mu)) +
209  (grad_x(lambda)+grad_u(lambda))*dgradu_dt(mu) );
210 
211  const libMesh::Real gamma_v = 0.5*( dgradv_dt(lambda)*(grad_y(mu)+grad_v(mu)) +
212  (grad_y(lambda)+grad_v(lambda))*dgradv_dt(mu) );
213 
214  const libMesh::Real x_term = C1*( (grad_x(beta)+grad_u(beta))*u_gradphi(alpha) +
215  (grad_x(alpha)+grad_u(alpha))*u_gradphi(beta) );
216 
217  const libMesh::Real y_term = C1*( (grad_y(beta)+grad_v(beta))*u_gradphi(alpha) +
218  (grad_y(alpha)+grad_v(alpha))*u_gradphi(beta) );
219 
220  libMesh::Real gamma_sum = gamma_u + gamma_v;
221  if( this->_disp_vars.dim() == 3 )
222  {
223  const libMesh::Real gamma_w = 0.5*( dgradw_dt(lambda)*(grad_z(mu)+grad_w(mu)) +
224  (grad_z(lambda)+grad_w(lambda))*dgradw_dt(mu) );
225  gamma_sum += gamma_w;
226  }
227 
228  Fu(i) += x_term*gamma_sum;
229 
230  Fv(i) += y_term*gamma_sum;
231 
232  if( this->_disp_vars.dim() == 3 )
233  {
234  const libMesh::Real z_term = C1*( (grad_z(beta)+grad_w(beta))*u_gradphi(alpha) +
235  (grad_z(alpha)+grad_w(alpha))*u_gradphi(beta) );
236 
237  (*Fw)(i) += z_term*gamma_sum;
238  }
239  }
240  }
241  }
242  }
243  }
244 
245  if( compute_jacobian )
246  {
247  for (unsigned int i=0; i != n_u_dofs; i++)
248  {
249  libMesh::RealGradient u_gradphi_i( dphi_dxi[i][qp], dphi_deta[i][qp] );
250 
251  for (unsigned int j=0; j != n_u_dofs; j++)
252  {
253  libMesh::RealGradient u_gradphi_j( dphi_dxi[j][qp], dphi_deta[j][qp] );
254 
255  for( unsigned int alpha = 0; alpha < manifold_dim; alpha++ )
256  {
257  for( unsigned int beta = 0; beta < manifold_dim; beta++ )
258  {
259  const libMesh::Real common_factor = _lambda_factor*0.5*this->_h0*jac;
260 
261  const libMesh::Real factor = common_factor*tau(alpha,beta);
262 
263  const libMesh::Real ddiagterm = u_gradphi_j(beta)*u_gradphi_i(alpha) + u_gradphi_j(alpha)*u_gradphi_i(beta);
264 
265  Kuu(i,j) += factor*ddiagterm*context.get_elem_solution_rate_derivative();
266  Kvv(i,j) += factor*ddiagterm*context.get_elem_solution_rate_derivative();
267 
268  if( this->_disp_vars.dim() == 3 )
269  (*Kww)(i,j) += factor*ddiagterm*context.get_elem_solution_rate_derivative();
270 
271  libMesh::Real u_diagterm = dgradu_dt(beta)*u_gradphi_i(alpha) + dgradu_dt(alpha)*u_gradphi_i(beta);
272  libMesh::Real v_diagterm = dgradv_dt(beta)*u_gradphi_i(alpha) + dgradv_dt(alpha)*u_gradphi_i(beta);
273  libMesh::Real w_diagterm = 0.0;
274  if( this->_disp_vars.dim() == 3 )
275  w_diagterm = dgradw_dt(beta)*u_gradphi_i(alpha) + dgradw_dt(alpha)*u_gradphi_i(beta);
276 
277  for( unsigned int lambda = 0; lambda < manifold_dim; lambda++ )
278  {
279  for( unsigned int mu = 0; mu < manifold_dim; mu++ )
280  {
281  const libMesh::Real dgamma_du = 0.5*( u_gradphi_j(lambda)*(grad_x(mu)+grad_u(mu)) +
282  (grad_x(lambda)+grad_u(lambda))*u_gradphi_j(mu) );
283 
284  const libMesh::Real dgamma_dv = 0.5*( u_gradphi_j(lambda)*(grad_y(mu)+grad_v(mu)) +
285  (grad_y(lambda)+grad_v(lambda))*u_gradphi_j(mu) );
286 
287  const libMesh::Real C1 = common_factor*C(alpha,beta,lambda,mu);
288 
289  const libMesh::Real dfactor_du = C1*dgamma_du*context.get_elem_solution_derivative();
290  const libMesh::Real dfactor_dv = C1*dgamma_dv*context.get_elem_solution_derivative();
291 
292  Kuu(i,j) += u_diagterm*dfactor_du;
293  Kuv(i,j) += u_diagterm*dfactor_dv;
294 
295  Kvu(i,j) += v_diagterm*dfactor_du;
296  Kvv(i,j) += v_diagterm*dfactor_dv;
297 
298  if( this->_disp_vars.dim() == 3 )
299  {
300  const libMesh::Real dgamma_dw = 0.5*( u_gradphi_j(lambda)*(grad_z(mu)+grad_w(mu)) +
301  (grad_z(lambda)+grad_w(lambda))*u_gradphi_j(mu) );
302 
303  const libMesh::Real dfactor_dw = C1*dgamma_dw*context.get_elem_solution_derivative();
304 
305  (*Kuw)(i,j) += u_diagterm*dfactor_dw;
306 
307  (*Kvw)(i,j) += v_diagterm*dfactor_dw;
308 
309  (*Kwu)(i,j) += w_diagterm*dfactor_du;
310  (*Kwv)(i,j) += w_diagterm*dfactor_dv;
311  (*Kww)(i,j) += w_diagterm*dfactor_dw;
312  }
313 
314  const libMesh::Real gamma_u = 0.5*( dgradu_dt(lambda)*(grad_x(mu)+grad_u(mu)) +
315  (grad_x(lambda)+grad_u(lambda))*dgradu_dt(mu) );
316 
317  const libMesh::Real gamma_v = 0.5*( dgradv_dt(lambda)*(grad_y(mu)+grad_v(mu)) +
318  (grad_y(lambda)+grad_v(lambda))*dgradv_dt(mu) );
319 
320  const libMesh::Real x_term = C1*( (grad_x(beta)+grad_u(beta))*u_gradphi_i(alpha) +
321  (grad_x(alpha)+grad_u(alpha))*u_gradphi_i(beta) );
322 
323  const libMesh::Real y_term = C1*( (grad_y(beta)+grad_v(beta))*u_gradphi_i(alpha) +
324  (grad_y(alpha)+grad_v(alpha))*u_gradphi_i(beta) );
325 
326  libMesh::Real gamma_sum = gamma_u + gamma_v;
327  if( this->_disp_vars.dim() == 3 )
328  {
329  const libMesh::Real gamma_w = 0.5*( dgradw_dt(lambda)*(grad_z(mu)+grad_w(mu)) +
330  (grad_z(lambda)+grad_w(lambda))*dgradw_dt(mu) );
331  gamma_sum += gamma_w;
332  }
333 
334  // Here, we're missing derivatives of C(alpha,beta,lambda,mu) w.r.t. strain
335  // Nonzero for hyperelasticity models
336  const libMesh::Real dxterm_du = C1*( u_gradphi_j(beta)*u_gradphi_i(alpha)
337  + u_gradphi_j(alpha)*u_gradphi_i(beta) )*context.get_elem_solution_derivative();
338 
339  const libMesh::Real dyterm_dv = dxterm_du;
340 
341  Kuu(i,j) += gamma_sum*dxterm_du;
342 
343  Kvv(i,j) += gamma_sum*dyterm_dv;
344 
345  libMesh::Real z_term = 0.0;
346 
347  if( this->_disp_vars.dim() == 3 )
348  {
349  z_term = C1*( (grad_z(beta)+grad_w(beta))*u_gradphi_i(alpha) +
350  (grad_z(alpha)+grad_w(alpha))*u_gradphi_i(beta) );
351 
352  // Here, we're missing derivatives of C(alpha,beta,lambda,mu) w.r.t. strain
353  // Nonzero for hyperelasticity models
354  const libMesh::Real dzterm_dw = dxterm_du;
355 
356  (*Kww)(i,j) += gamma_sum*dzterm_dw;
357  }
358 
359  const libMesh::Real dgamma_sum_du = 0.5*( u_gradphi_j(lambda)*(grad_x(mu)+grad_u(mu))*context.get_elem_solution_rate_derivative()
360  + dgradu_dt(lambda)*u_gradphi_j(mu)*context.get_elem_solution_derivative()
361  + u_gradphi_j(mu)*(grad_x(lambda)+grad_u(lambda))*context.get_elem_solution_rate_derivative()
362  + dgradu_dt(mu)*u_gradphi_j(lambda)*context.get_elem_solution_derivative() );
363 
364  const libMesh::Real dgamma_sum_dv = 0.5*( u_gradphi_j(lambda)*(grad_y(mu)+grad_v(mu))*context.get_elem_solution_rate_derivative()
365  + dgradv_dt(lambda)*u_gradphi_j(mu)*context.get_elem_solution_derivative()
366  + u_gradphi_j(mu)*(grad_y(lambda)+grad_v(lambda))*context.get_elem_solution_rate_derivative()
367  + dgradv_dt(mu)*u_gradphi_j(lambda)*context.get_elem_solution_derivative() );
368 
369  Kuu(i,j) += x_term*dgamma_sum_du;
370  Kuv(i,j) += x_term*dgamma_sum_dv;
371 
372  Kvu(i,j) += y_term*dgamma_sum_du;
373  Kvv(i,j) += y_term*dgamma_sum_dv;
374 
375  if( this->_disp_vars.dim() == 3 )
376  {
377  const libMesh::Real dgamma_sum_dw = 0.5*( u_gradphi_j(lambda)*(grad_z(mu)+grad_w(mu))*context.get_elem_solution_rate_derivative()
378  + dgradw_dt(lambda)*u_gradphi_j(mu)*context.get_elem_solution_derivative()
379  + u_gradphi_j(mu)*(grad_z(lambda)+grad_w(lambda))*context.get_elem_solution_rate_derivative()
380  + dgradw_dt(mu)*u_gradphi_j(lambda)*context.get_elem_solution_derivative() );
381 
382  (*Kuw)(i,j) += x_term*dgamma_sum_dw;
383  (*Kvw)(i,j) += y_term*dgamma_sum_dw;
384  (*Kwu)(i,j) += z_term*dgamma_sum_du;
385  (*Kwv)(i,j) += z_term*dgamma_sum_dv;
386  (*Kww)(i,j) += z_term*dgamma_sum_dw;
387  }
388  }
389  }
390  }
391  }
392  }
393  }
394  }
395 
396  }
397  }
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.
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)

Member Data Documentation

template<typename StressStrainLaw >
libMesh::Real GRINS::ElasticMembraneRayleighDamping< StressStrainLaw >::_lambda_factor
protected

Definition at line 49 of file elastic_membrane_rayleigh_damping.h.

template<typename StressStrainLaw >
libMesh::Real GRINS::ElasticMembraneRayleighDamping< StressStrainLaw >::_mu_factor
protected

Definition at line 50 of file elastic_membrane_rayleigh_damping.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