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

#include <elastic_cable_rayleigh_damping.h>

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

Public Member Functions

 ElasticCableRayleighDamping (const PhysicsName &physics_name, const GetPot &input, bool is_compressible)
 
virtual ~ElasticCableRayleighDamping ()
 
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::ElasticCableBase< StressStrainLaw >
 ElasticCableBase (const GRINS::PhysicsName &physics_name, const GetPot &input, bool is_compressible)
 
virtual ~ElasticCableBase ()
 
- Public Member Functions inherited from GRINS::ElasticCableAbstract
 ElasticCableAbstract (const GRINS::PhysicsName &physics_name, const GetPot &input)
 
virtual ~ElasticCableAbstract ()
 
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 void init_variables (libMesh::FEMSystem *)
 Initialize variables for this physics. 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 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::ElasticCableBase< StressStrainLaw >
StressStrainLaw _stress_strain_law
 
bool _is_compressible
 
- Protected Attributes inherited from GRINS::ElasticCableAbstract
libMesh::Real _A
 Cross-sectional area of the cable. More...
 
libMesh::Real _rho
 Cable density. 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

 ElasticCableRayleighDamping ()
 

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::ElasticCableBase< 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::ElasticCableAbstract
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::ElasticCableRayleighDamping< StressStrainLaw >

Definition at line 33 of file elastic_cable_rayleigh_damping.h.

Constructor & Destructor Documentation

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

Definition at line 39 of file elastic_cable_rayleigh_damping.C.

References GRINS::PhysicsNaming::elastic_cable(), GRINS::PhysicsNaming::elastic_cable_rayleigh_damping(), and GRINS::Physics::parse_enabled_subdomains().

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

Definition at line 41 of file elastic_cable_rayleigh_damping.h.

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

Member Function Documentation

template<typename StressStrainLaw >
void GRINS::ElasticCableRayleighDamping< 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_cable_rayleigh_damping.C.

67  {
68  // First, do the "mass" contribution
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 contribution
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, Jacobians that we're populating
82  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_disp_vars.u());
83  libMesh::DenseSubVector<libMesh::Number>* Fv = NULL;
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 = NULL;
88  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
89 
90  libMesh::DenseSubMatrix<libMesh::Number>* Kvu = NULL;
91  libMesh::DenseSubMatrix<libMesh::Number>* Kvv = NULL;
92  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
93 
94  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
95  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
96  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
97 
98  if( this->_disp_vars.dim() >= 2 )
99  {
100  Fv = &context.get_elem_residual(this->_disp_vars.v());
101 
102  Kuv = &context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.v());
103  Kvu = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.u());
104  Kvv = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.v());
105  }
106 
107  if( this->_disp_vars.dim() == 3 )
108  {
109  Fw = &context.get_elem_residual(this->_disp_vars.w());
110 
111  Kuw = &context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.w());
112  Kvw = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.w());
113  Kwu = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.u());
114  Kwv = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.v());
115  Kww = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.w());
116  }
117 
118  unsigned int n_qpoints = context.get_element_qrule().n_points();
119 
120  // All shape function gradients are w.r.t. master element coordinates
121  const std::vector<std::vector<libMesh::Real> >& dphi_dxi = this->get_fe(context)->get_dphidxi();
122 
123  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
124  const libMesh::DenseSubVector<libMesh::Number>* v_coeffs = NULL;
125  const libMesh::DenseSubVector<libMesh::Number>* w_coeffs = NULL;
126 
127  if( this->_disp_vars.dim() >= 2 )
128  v_coeffs = &context.get_elem_solution( this->_disp_vars.v() );
129 
130  if( this->_disp_vars.dim() == 3 )
131  w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
132 
133  const libMesh::DenseSubVector<libMesh::Number>& dudt_coeffs =
134  context.get_elem_solution_rate( this->_disp_vars.u() );
135  const libMesh::DenseSubVector<libMesh::Number>* dvdt_coeffs = NULL;
136  const libMesh::DenseSubVector<libMesh::Number>* dwdt_coeffs = NULL;
137 
138  if( this->_disp_vars.dim() >= 2 )
139  dvdt_coeffs = &context.get_elem_solution_rate( this->_disp_vars.v() );
140 
141  if( this->_disp_vars.dim() == 3 )
142  dwdt_coeffs = &context.get_elem_solution_rate( this->_disp_vars.w() );
143 
144  // Need these to build up the covariant and contravariant metric tensors
145  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
146 
147  const unsigned int dim = 1; // The cable dimension is always 1 for this physics
148 
149  for (unsigned int qp=0; qp != n_qpoints; qp++)
150  {
151  // Gradients are w.r.t. master element coordinates
152  libMesh::Gradient grad_u, grad_v, grad_w;
153  libMesh::Gradient dgradu_dt, dgradv_dt, dgradw_dt;
154 
155  for( unsigned int d = 0; d < n_u_dofs; d++ )
156  {
157  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp] );
158  grad_u += u_coeffs(d)*u_gradphi;
159  dgradu_dt += dudt_coeffs(d)*u_gradphi;
160 
161  if( this->_disp_vars.dim() >= 2 )
162  {
163  grad_v += (*v_coeffs)(d)*u_gradphi;
164  dgradv_dt += (*dvdt_coeffs)(d)*u_gradphi;
165  }
166 
167  if( this->_disp_vars.dim() == 3 )
168  {
169  grad_w += (*w_coeffs)(d)*u_gradphi;
170  dgradw_dt += (*dwdt_coeffs)(d)*u_gradphi;
171  }
172  }
173 
174  libMesh::RealGradient grad_x( dxdxi[qp](0) );
175  libMesh::RealGradient grad_y( dxdxi[qp](1) );
176  libMesh::RealGradient grad_z( dxdxi[qp](2) );
177 
178  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
179  libMesh::Real lambda_sq = 0;
180 
181  this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
182  grad_u, grad_v, grad_w,
183  a_cov, a_contra, A_cov, A_contra,
184  lambda_sq );
185 
186  // Compute stress tensor
188  ElasticityTensor C;
189  this->_stress_strain_law.compute_stress_and_elasticity(dim,a_contra,a_cov,A_contra,A_cov,tau,C);
190 
191  libMesh::Real jac = JxW[qp];
192 
193  for (unsigned int i=0; i != n_u_dofs; i++)
194  {
195  libMesh::RealGradient u_gradphi( dphi_dxi[i][qp] );
196 
197  libMesh::Real u_diag_factor = _lambda_factor*this->_A*jac*tau(0,0)*dgradu_dt(0)*u_gradphi(0);
198  libMesh::Real v_diag_factor = _lambda_factor*this->_A*jac*tau(0,0)*dgradv_dt(0)*u_gradphi(0);
199  libMesh::Real w_diag_factor = _lambda_factor*this->_A*jac*tau(0,0)*dgradw_dt(0)*u_gradphi(0);
200 
201  const libMesh::Real C1 = _lambda_factor*this->_A*jac*C(0,0,0,0)*u_gradphi(0);
202 
203  const libMesh::Real gamma_u = (grad_x(0)+grad_u(0));
204 
205  libMesh::Real gamma_v = 0.0;
206  libMesh::Real gamma_w = 0.0;
207 
208  if( this->_disp_vars.dim() >= 2 )
209  gamma_v = (grad_y(0)+grad_v(0));
210 
211  if( this->_disp_vars.dim() == 3 )
212  gamma_w = (grad_z(0)+grad_w(0));
213 
214  const libMesh::Real x_term = C1*gamma_u;
215  const libMesh::Real y_term = C1*gamma_v;
216  const libMesh::Real z_term = C1*gamma_w;
217 
218  libMesh::Real dt_term = dgradu_dt(0)*gamma_u + dgradv_dt(0)*gamma_v + dgradw_dt(0)*gamma_w;
219 
220  Fu(i) += u_diag_factor + x_term*dt_term;
221 
222  if( this->_disp_vars.dim() >= 2 )
223  (*Fv)(i) += v_diag_factor + y_term*dt_term;
224 
225  if( this->_disp_vars.dim() == 3 )
226  (*Fw)(i) += w_diag_factor + z_term*dt_term;
227  }
228 
229  if( compute_jacobian )
230  {
231  for(unsigned int i=0; i != n_u_dofs; i++)
232  {
233  libMesh::RealGradient u_gradphi_I( dphi_dxi[i][qp] );
234 
235  for(unsigned int j=0; j != n_u_dofs; j++)
236  {
237  libMesh::RealGradient u_gradphi_J( dphi_dxi[j][qp] );
238 
239  libMesh::Real common_factor = _lambda_factor*this->_A*jac*u_gradphi_I(0);
240 
241  const libMesh::Real diag_term_1 = common_factor*tau(0,0)*u_gradphi_J(0)*context.get_elem_solution_rate_derivative();
242 
243  const libMesh::Real dgamma_du = ( u_gradphi_J(0)*(grad_x(0)+grad_u(0)) );
244 
245  const libMesh::Real dgamma_dv = ( u_gradphi_J(0)*(grad_y(0)+grad_v(0)) );
246 
247  const libMesh::Real dgamma_dw = ( u_gradphi_J(0)*(grad_z(0)+grad_w(0)) );
248 
249  const libMesh::Real diag_term_2_factor = common_factor*C(0,0,0,0)*context.get_elem_solution_derivative();
250 
251  Kuu(i,j) += diag_term_1 + dgradu_dt(0)*diag_term_2_factor*dgamma_du;
252 
253  if( this->_disp_vars.dim() >= 2 )
254  {
255  (*Kuv)(i,j) += dgradu_dt(0)*diag_term_2_factor*dgamma_dv;
256  (*Kvu)(i,j) += dgradv_dt(0)*diag_term_2_factor*dgamma_du;
257  (*Kvv)(i,j) += diag_term_1 + dgradv_dt(0)*diag_term_2_factor*dgamma_dv;
258  }
259 
260  if( this->_disp_vars.dim() == 3 )
261  {
262  (*Kuw)(i,j) += dgradu_dt(0)*diag_term_2_factor*dgamma_dw;
263  (*Kvw)(i,j) += dgradv_dt(0)*diag_term_2_factor*dgamma_dw;
264  (*Kwu)(i,j) += dgradw_dt(0)*diag_term_2_factor*dgamma_du;
265  (*Kwv)(i,j) += dgradw_dt(0)*diag_term_2_factor*dgamma_dv;
266  (*Kww)(i,j) += diag_term_1 + dgradw_dt(0)*diag_term_2_factor*dgamma_dw;
267  }
268 
269  const libMesh::Real C1 = common_factor*C(0,0,0,0);
270 
271  const libMesh::Real gamma_u = (grad_x(0)+grad_u(0));
272  const libMesh::Real gamma_v = (grad_y(0)+grad_v(0));
273  const libMesh::Real gamma_w = (grad_z(0)+grad_w(0));
274 
275  const libMesh::Real x_term = C1*gamma_u;
276  const libMesh::Real y_term = C1*gamma_v;
277  const libMesh::Real z_term = C1*gamma_w;
278 
279  const libMesh::Real ddtterm_du = u_gradphi_J(0)*(gamma_u*context.get_elem_solution_rate_derivative()
280  + dgradu_dt(0)*context.get_elem_solution_derivative());
281 
282  libMesh::Real ddtterm_dv = 0.0;
283  if( this->_disp_vars.dim() >= 2 )
284  ddtterm_dv = u_gradphi_J(0)*(gamma_v*context.get_elem_solution_rate_derivative()
285  + dgradv_dt(0)*context.get_elem_solution_derivative());
286 
287  libMesh::Real ddtterm_dw = 0.0;
288  if( this->_disp_vars.dim() == 3 )
289  ddtterm_dw = u_gradphi_J(0)*(gamma_w*context.get_elem_solution_rate_derivative()
290  + dgradw_dt(0)*context.get_elem_solution_derivative());
291 
292  Kuu(i,j) += x_term*ddtterm_du;
293 
294  if( this->_disp_vars.dim() >= 2 )
295  {
296  (*Kuv)(i,j) += x_term*ddtterm_dv;
297  (*Kvu)(i,j) += y_term*ddtterm_du;
298  (*Kvv)(i,j) += y_term*ddtterm_dv;
299  }
300 
301  if( this->_disp_vars.dim() == 3 )
302  {
303  (*Kuw)(i,j) += x_term*ddtterm_dw;
304  (*Kvw)(i,j) += y_term*ddtterm_dw;
305  (*Kwu)(i,j) += z_term*ddtterm_du;
306  (*Kwv)(i,j) += z_term*ddtterm_dv;
307  (*Kww)(i,j) += z_term*ddtterm_dw;
308  }
309 
310  libMesh::Real dt_term = dgradu_dt(0)*gamma_u + dgradv_dt(0)*gamma_v + dgradw_dt(0)*gamma_w;
311 
312  // Here, we're missing derivatives of C(0,0,0,0) w.r.t. strain
313  // Nonzero for hyperelasticity models
314  const libMesh::Real dxterm_du = C1*u_gradphi_J(0)*context.get_elem_solution_derivative();
315  const libMesh::Real dyterm_dv = dxterm_du;
316  const libMesh::Real dzterm_dw = dxterm_du;
317 
318  Kuu(i,j) += dxterm_du*dt_term;
319 
320  if( this->_disp_vars.dim() >= 2 )
321  (*Kvv)(i,j) += dyterm_dv*dt_term;
322 
323  if( this->_disp_vars.dim() == 3 )
324  (*Kww)(i,j) += dzterm_dw*dt_term;
325 
326  } // end j-loop
327  } // end i-loop
328  } // end if(compute_jacobian)
329  } // end qp loop
330  }
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)
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 mass_residual_impl(bool compute_jacobian, AssemblyContext &context, InteriorFuncType interior_solution, VarDerivType get_solution_deriv, libMesh::Real mu=1.0)
Implementation of mass_residual.
StressStrainLaw _stress_strain_law
libMesh::Real _A
Cross-sectional area of the cable.
unsigned int dim() const
Number of components.
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)

Member Data Documentation

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

Definition at line 49 of file elastic_cable_rayleigh_damping.h.

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

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