33 #include "libmesh/getpot.h"
34 #include "libmesh/quadrature.h"
38 template<
typename StressStrainLaw>
41 bool 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))
56 libmesh_error_msg(
"ERROR: Cannot specify subdomains for "
58 +
"! Must specify subdomains through "
64 template<
typename StressStrainLaw>
69 this->mass_residual_impl(compute_jacobian,
71 &libMesh::FEMContext::interior_rate,
72 &libMesh::DiffContext::get_elem_solution_rate_derivative,
76 const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
78 const std::vector<libMesh::Real> &JxW =
79 this->get_fe(context)->get_JxW();
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;
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;
90 libMesh::DenseSubMatrix<libMesh::Number>* Kvu = NULL;
91 libMesh::DenseSubMatrix<libMesh::Number>* Kvv = NULL;
92 libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
94 libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
95 libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
96 libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
98 if( this->_disp_vars.dim() >= 2 )
100 Fv = &context.get_elem_residual(this->_disp_vars.v());
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());
107 if( this->_disp_vars.dim() == 3 )
109 Fw = &context.get_elem_residual(this->_disp_vars.w());
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());
118 unsigned int n_qpoints = context.get_element_qrule().n_points();
121 const std::vector<std::vector<libMesh::Real> >& dphi_dxi = this->get_fe(context)->get_dphidxi();
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;
127 if( this->_disp_vars.dim() >= 2 )
128 v_coeffs = &context.get_elem_solution( this->_disp_vars.v() );
130 if( this->_disp_vars.dim() == 3 )
131 w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
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;
138 if( this->_disp_vars.dim() >= 2 )
139 dvdt_coeffs = &context.get_elem_solution_rate( this->_disp_vars.v() );
141 if( this->_disp_vars.dim() == 3 )
142 dwdt_coeffs = &context.get_elem_solution_rate( this->_disp_vars.w() );
145 const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
147 const unsigned int dim = 1;
149 for (
unsigned int qp=0; qp != n_qpoints; qp++)
152 libMesh::Gradient grad_u, grad_v, grad_w;
153 libMesh::Gradient dgradu_dt, dgradv_dt, dgradw_dt;
155 for(
unsigned int d = 0; d < n_u_dofs; d++ )
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;
161 if( this->_disp_vars.dim() >= 2 )
163 grad_v += (*v_coeffs)(d)*u_gradphi;
164 dgradv_dt += (*dvdt_coeffs)(d)*u_gradphi;
167 if( this->_disp_vars.dim() == 3 )
169 grad_w += (*w_coeffs)(d)*u_gradphi;
170 dgradw_dt += (*dwdt_coeffs)(d)*u_gradphi;
174 libMesh::RealGradient grad_x( dxdxi[qp](0) );
175 libMesh::RealGradient grad_y( dxdxi[qp](1) );
176 libMesh::RealGradient grad_z( dxdxi[qp](2) );
179 libMesh::Real lambda_sq = 0;
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,
189 this->_stress_strain_law.compute_stress_and_elasticity(dim,a_contra,a_cov,A_contra,A_cov,tau,C);
191 libMesh::Real jac = JxW[qp];
193 for (
unsigned int i=0; i != n_u_dofs; i++)
195 libMesh::RealGradient u_gradphi( dphi_dxi[i][qp] );
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);
201 const libMesh::Real C1 = _lambda_factor*this->_A*jac*C(0,0,0,0)*u_gradphi(0);
203 const libMesh::Real gamma_u = (grad_x(0)+grad_u(0));
205 libMesh::Real gamma_v = 0.0;
206 libMesh::Real gamma_w = 0.0;
208 if( this->_disp_vars.dim() >= 2 )
209 gamma_v = (grad_y(0)+grad_v(0));
211 if( this->_disp_vars.dim() == 3 )
212 gamma_w = (grad_z(0)+grad_w(0));
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;
218 libMesh::Real dt_term = dgradu_dt(0)*gamma_u + dgradv_dt(0)*gamma_v + dgradw_dt(0)*gamma_w;
220 Fu(i) += u_diag_factor + x_term*dt_term;
222 if( this->_disp_vars.dim() >= 2 )
223 (*Fv)(i) += v_diag_factor + y_term*dt_term;
225 if( this->_disp_vars.dim() == 3 )
226 (*Fw)(i) += w_diag_factor + z_term*dt_term;
229 if( compute_jacobian )
231 for(
unsigned int i=0; i != n_u_dofs; i++)
233 libMesh::RealGradient u_gradphi_I( dphi_dxi[i][qp] );
235 for(
unsigned int j=0; j != n_u_dofs; j++)
237 libMesh::RealGradient u_gradphi_J( dphi_dxi[j][qp] );
239 libMesh::Real common_factor = _lambda_factor*this->_A*jac*u_gradphi_I(0);
241 const libMesh::Real diag_term_1 = common_factor*tau(0,0)*u_gradphi_J(0)*context.get_elem_solution_rate_derivative();
243 const libMesh::Real dgamma_du = ( u_gradphi_J(0)*(grad_x(0)+grad_u(0)) );
245 const libMesh::Real dgamma_dv = ( u_gradphi_J(0)*(grad_y(0)+grad_v(0)) );
247 const libMesh::Real dgamma_dw = ( u_gradphi_J(0)*(grad_z(0)+grad_w(0)) );
249 const libMesh::Real diag_term_2_factor = common_factor*C(0,0,0,0)*context.get_elem_solution_derivative();
251 Kuu(i,j) += diag_term_1 + dgradu_dt(0)*diag_term_2_factor*dgamma_du;
253 if( this->_disp_vars.dim() >= 2 )
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;
260 if( this->_disp_vars.dim() == 3 )
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;
269 const libMesh::Real C1 = common_factor*C(0,0,0,0);
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));
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;
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());
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());
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());
292 Kuu(i,j) += x_term*ddtterm_du;
294 if( this->_disp_vars.dim() >= 2 )
296 (*Kuv)(i,j) += x_term*ddtterm_dv;
297 (*Kvu)(i,j) += y_term*ddtterm_du;
298 (*Kvv)(i,j) += y_term*ddtterm_dv;
301 if( this->_disp_vars.dim() == 3 )
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;
310 libMesh::Real dt_term = dgradu_dt(0)*gamma_u + dgradv_dt(0)*gamma_v + dgradw_dt(0)*gamma_w;
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;
318 Kuu(i,j) += dxterm_du*dt_term;
320 if( this->_disp_vars.dim() >= 2 )
321 (*Kvv)(i,j) += dyterm_dv*dt_term;
323 if( this->_disp_vars.dim() == 3 )
324 (*Kww)(i,j) += dzterm_dw*dt_term;
static PhysicsName elastic_cable_rayleigh_damping()
static PhysicsName elastic_cable()
virtual void damping_residual(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
void parse_enabled_subdomains(const GetPot &input, const std::string &physics_name)
ElasticCableRayleighDamping()