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>
70 this->mass_residual_impl(compute_jacobian,
72 &libMesh::FEMContext::interior_rate,
73 &libMesh::DiffContext::get_elem_solution_rate_derivative,
77 const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
79 const std::vector<libMesh::Real> &JxW =
80 this->get_fe(context)->get_JxW();
83 libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_disp_vars.u());
84 libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(this->_disp_vars.v());
85 libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(this->_disp_vars.w());
89 libMesh::DenseSubMatrix<libMesh::Number> &Kuu = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.u());
90 libMesh::DenseSubMatrix<libMesh::Number> &Kuv = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.v());
91 libMesh::DenseSubMatrix<libMesh::Number> &Kuw = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.w());
92 libMesh::DenseSubMatrix<libMesh::Number> &Kvu = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.u());
93 libMesh::DenseSubMatrix<libMesh::Number> &Kvv = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.v());
94 libMesh::DenseSubMatrix<libMesh::Number> &Kvw = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.w());
95 libMesh::DenseSubMatrix<libMesh::Number> &Kwu = context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.u());
96 libMesh::DenseSubMatrix<libMesh::Number> &Kwv = context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.v());
97 libMesh::DenseSubMatrix<libMesh::Number> &Kww = context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.w());
99 unsigned int n_qpoints = context.get_element_qrule().n_points();
102 const std::vector<std::vector<libMesh::Real> >& dphi_dxi = this->get_fe(context)->get_dphidxi();
104 const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
105 const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( this->_disp_vars.v() );
106 const libMesh::DenseSubVector<libMesh::Number>& w_coeffs = context.get_elem_solution( this->_disp_vars.w() );
108 const libMesh::DenseSubVector<libMesh::Number>& dudt_coeffs = context.get_elem_solution_rate( this->_disp_vars.u() );
109 const libMesh::DenseSubVector<libMesh::Number>& dvdt_coeffs = context.get_elem_solution_rate( this->_disp_vars.v() );
110 const libMesh::DenseSubVector<libMesh::Number>& dwdt_coeffs = context.get_elem_solution_rate( this->_disp_vars.w() );
113 const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
115 const unsigned int dim = 1;
117 for (
unsigned int qp=0; qp != n_qpoints; qp++)
120 libMesh::Gradient grad_u, grad_v, grad_w;
121 libMesh::Gradient dgradu_dt, dgradv_dt, dgradw_dt;
123 for(
unsigned int d = 0; d < n_u_dofs; d++ )
125 libMesh::RealGradient u_gradphi( dphi_dxi[d][qp] );
126 grad_u += u_coeffs(d)*u_gradphi;
127 grad_v += v_coeffs(d)*u_gradphi;
128 grad_w += w_coeffs(d)*u_gradphi;
130 dgradu_dt += dudt_coeffs(d)*u_gradphi;
131 dgradv_dt += dvdt_coeffs(d)*u_gradphi;
132 dgradw_dt += dwdt_coeffs(d)*u_gradphi;
135 libMesh::RealGradient grad_x( dxdxi[qp](0) );
136 libMesh::RealGradient grad_y( dxdxi[qp](1) );
137 libMesh::RealGradient grad_z( dxdxi[qp](2) );
140 libMesh::Real lambda_sq = 0;
142 this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
143 grad_u, grad_v, grad_w,
144 a_cov, a_contra, A_cov, A_contra,
150 this->_stress_strain_law.compute_stress_and_elasticity(dim,a_contra,a_cov,A_contra,A_cov,tau,C);
152 libMesh::Real jac = JxW[qp];
154 for (
unsigned int i=0; i != n_u_dofs; i++)
156 libMesh::RealGradient u_gradphi( dphi_dxi[i][qp] );
158 libMesh::Real u_diag_factor = _lambda_factor*this->_A*jac*tau(0,0)*dgradu_dt(0)*u_gradphi(0);
159 libMesh::Real v_diag_factor = _lambda_factor*this->_A*jac*tau(0,0)*dgradv_dt(0)*u_gradphi(0);
160 libMesh::Real w_diag_factor = _lambda_factor*this->_A*jac*tau(0,0)*dgradw_dt(0)*u_gradphi(0);
162 const libMesh::Real C1 = _lambda_factor*this->_A*jac*C(0,0,0,0)*u_gradphi(0);
164 const libMesh::Real gamma_u = (grad_x(0)+grad_u(0));
165 const libMesh::Real gamma_v = (grad_y(0)+grad_v(0));
166 const libMesh::Real gamma_w = (grad_z(0)+grad_w(0));
168 const libMesh::Real x_term = C1*gamma_u;
169 const libMesh::Real y_term = C1*gamma_v;
170 const libMesh::Real z_term = C1*gamma_w;
172 const libMesh::Real dt_term = dgradu_dt(0)*gamma_u + dgradv_dt(0)*gamma_v + dgradw_dt(0)*gamma_w;
174 Fu(i) += u_diag_factor + x_term*dt_term;
175 Fv(i) += v_diag_factor + y_term*dt_term;
176 Fw(i) += w_diag_factor + z_term*dt_term;
179 if( compute_jacobian )
181 for(
unsigned int i=0; i != n_u_dofs; i++)
183 libMesh::RealGradient u_gradphi_I( dphi_dxi[i][qp] );
185 for(
unsigned int j=0; j != n_u_dofs; j++)
187 libMesh::RealGradient u_gradphi_J( dphi_dxi[j][qp] );
189 libMesh::Real common_factor = _lambda_factor*this->_A*jac*u_gradphi_I(0);
191 const libMesh::Real diag_term_1 = common_factor*tau(0,0)*u_gradphi_J(0)*context.get_elem_solution_rate_derivative();
193 const libMesh::Real dgamma_du = ( u_gradphi_J(0)*(grad_x(0)+grad_u(0)) );
195 const libMesh::Real dgamma_dv = ( u_gradphi_J(0)*(grad_y(0)+grad_v(0)) );
197 const libMesh::Real dgamma_dw = ( u_gradphi_J(0)*(grad_z(0)+grad_w(0)) );
199 const libMesh::Real diag_term_2_factor = common_factor*C(0,0,0,0)*context.get_elem_solution_derivative();
201 Kuu(i,j) += diag_term_1 + dgradu_dt(0)*diag_term_2_factor*dgamma_du;
202 Kuv(i,j) += dgradu_dt(0)*diag_term_2_factor*dgamma_dv;
203 Kuw(i,j) += dgradu_dt(0)*diag_term_2_factor*dgamma_dw;
205 Kvu(i,j) += dgradv_dt(0)*diag_term_2_factor*dgamma_du;
206 Kvv(i,j) += diag_term_1 + dgradv_dt(0)*diag_term_2_factor*dgamma_dv;
207 Kvw(i,j) += dgradv_dt(0)*diag_term_2_factor*dgamma_dw;
209 Kwu(i,j) += dgradw_dt(0)*diag_term_2_factor*dgamma_du;
210 Kwv(i,j) += dgradw_dt(0)*diag_term_2_factor*dgamma_dv;
211 Kww(i,j) += diag_term_1 + dgradw_dt(0)*diag_term_2_factor*dgamma_dw;
213 const libMesh::Real C1 = common_factor*C(0,0,0,0);
215 const libMesh::Real gamma_u = (grad_x(0)+grad_u(0));
216 const libMesh::Real gamma_v = (grad_y(0)+grad_v(0));
217 const libMesh::Real gamma_w = (grad_z(0)+grad_w(0));
219 const libMesh::Real x_term = C1*gamma_u;
220 const libMesh::Real y_term = C1*gamma_v;
221 const libMesh::Real z_term = C1*gamma_w;
223 const libMesh::Real ddtterm_du = u_gradphi_J(0)*(gamma_u*context.get_elem_solution_rate_derivative()
224 + dgradu_dt(0)*context.get_elem_solution_derivative());
226 const libMesh::Real ddtterm_dv = u_gradphi_J(0)*(gamma_v*context.get_elem_solution_rate_derivative()
227 + dgradv_dt(0)*context.get_elem_solution_derivative());
229 const libMesh::Real ddtterm_dw = u_gradphi_J(0)*(gamma_w*context.get_elem_solution_rate_derivative()
230 + dgradw_dt(0)*context.get_elem_solution_derivative());
232 Kuu(i,j) += x_term*ddtterm_du;
233 Kuv(i,j) += x_term*ddtterm_dv;
234 Kuw(i,j) += x_term*ddtterm_dw;
236 Kvu(i,j) += y_term*ddtterm_du;
237 Kvv(i,j) += y_term*ddtterm_dv;
238 Kvw(i,j) += y_term*ddtterm_dw;
240 Kwu(i,j) += z_term*ddtterm_du;
241 Kwv(i,j) += z_term*ddtterm_dv;
242 Kww(i,j) += z_term*ddtterm_dw;
244 const libMesh::Real dt_term = dgradu_dt(0)*gamma_u + dgradv_dt(0)*gamma_v + dgradw_dt(0)*gamma_w;
248 const libMesh::Real dxterm_du = C1*u_gradphi_J(0)*context.get_elem_solution_derivative();
249 const libMesh::Real dyterm_dv = dxterm_du;
250 const libMesh::Real dzterm_dw = dxterm_du;
252 Kuu(i,j) += dxterm_du*dt_term;
253 Kvv(i,j) += dyterm_dv*dt_term;
254 Kww(i,j) += dzterm_dw*dt_term;
static PhysicsName elastic_cable_rayleigh_damping()
static PhysicsName elastic_cable()
void parse_enabled_subdomains(const GetPot &input, const std::string &physics_name)
virtual void damping_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &)
Time dependent part(s) of physics for element interiors.
ElasticCableRayleighDamping()