GRINS-0.7.0
elastic_membrane_rayleigh_damping.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2016 Paul T. Bauman, Roy H. Stogner
7 // Copyright (C) 2010-2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 // This class
27 
28 // GRINS
29 #include "grins/physics_naming.h"
31 #include "grins/assembly_context.h"
32 
33 // libMesh
34 #include "libmesh/getpot.h"
35 #include "libmesh/quadrature.h"
36 
37 namespace GRINS
38 {
39  template<typename StressStrainLaw>
41  const GetPot& input,
42  bool is_compressible )
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  }
63 
64  template<typename StressStrainLaw>
66  AssemblyContext& context,
67  CachedValues& /*cache*/ )
68  {
69  // First do the mass part
70  this->mass_residual_impl(compute_jacobian,
71  context,
72  &libMesh::FEMContext::interior_rate,
73  &libMesh::DiffContext::get_elem_solution_rate_derivative,
74  _mu_factor);
75 
76  // Now do the stiffness part
77  const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
78 
79  const std::vector<libMesh::Real> &JxW =
80  this->get_fe(context)->get_JxW();
81 
82  // Residuals that we're populating
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());
86 
87  libMesh::DenseSubMatrix<libMesh::Number>& Kuu = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.u());
88  libMesh::DenseSubMatrix<libMesh::Number>& Kuv = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.v());
89  libMesh::DenseSubMatrix<libMesh::Number>& Kuw = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.w());
90 
91  libMesh::DenseSubMatrix<libMesh::Number>& Kvu = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.u());
92  libMesh::DenseSubMatrix<libMesh::Number>& Kvv = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.v());
93  libMesh::DenseSubMatrix<libMesh::Number>& Kvw = context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.w());
94 
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());
98 
99  unsigned int n_qpoints = context.get_element_qrule().n_points();
100 
101  // All shape function gradients are w.r.t. master element coordinates
102  const std::vector<std::vector<libMesh::Real> >& dphi_dxi =
103  this->get_fe(context)->get_dphidxi();
104 
105  const std::vector<std::vector<libMesh::Real> >& dphi_deta =
106  this->get_fe(context)->get_dphideta();
107 
108  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
109  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( this->_disp_vars.v() );
110  const libMesh::DenseSubVector<libMesh::Number>& w_coeffs = context.get_elem_solution( this->_disp_vars.w() );
111 
112  const libMesh::DenseSubVector<libMesh::Number>& dudt_coeffs = context.get_elem_solution_rate( this->_disp_vars.u() );
113  const libMesh::DenseSubVector<libMesh::Number>& dvdt_coeffs = context.get_elem_solution_rate( this->_disp_vars.v() );
114  const libMesh::DenseSubVector<libMesh::Number>& dwdt_coeffs = context.get_elem_solution_rate( this->_disp_vars.w() );
115 
116  // Need these to build up the covariant and contravariant metric tensors
117  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
118  const std::vector<libMesh::RealGradient>& dxdeta = this->get_fe(context)->get_dxyzdeta();
119 
120  const unsigned int dim = 2; // The manifold dimension is always 2 for this physics
121 
122  for (unsigned int qp=0; qp != n_qpoints; qp++)
123  {
124  // Gradients are w.r.t. master element coordinates
125  libMesh::Gradient grad_u, grad_v, grad_w;
126  libMesh::Gradient dgradu_dt, dgradv_dt, dgradw_dt;
127 
128  for( unsigned int d = 0; d < n_u_dofs; d++ )
129  {
130  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp], dphi_deta[d][qp] );
131  grad_u += u_coeffs(d)*u_gradphi;
132  grad_v += v_coeffs(d)*u_gradphi;
133  grad_w += w_coeffs(d)*u_gradphi;
134 
135  dgradu_dt += dudt_coeffs(d)*u_gradphi;
136  dgradv_dt += dvdt_coeffs(d)*u_gradphi;
137  dgradw_dt += dwdt_coeffs(d)*u_gradphi;
138  }
139 
140  libMesh::RealGradient grad_x( dxdxi[qp](0), dxdeta[qp](0) );
141  libMesh::RealGradient grad_y( dxdxi[qp](1), dxdeta[qp](1) );
142  libMesh::RealGradient grad_z( dxdxi[qp](2), dxdeta[qp](2) );
143 
144 
145  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
146  libMesh::Real lambda_sq = 0;
147 
148  this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
149  grad_u, grad_v, grad_w,
150  a_cov, a_contra, A_cov, A_contra,
151  lambda_sq );
152 
153  // Compute stress and elasticity tensors
156  this->_stress_strain_law.compute_stress_and_elasticity(dim,a_contra,a_cov,A_contra,A_cov,tau,C);
157 
158  libMesh::Real jac = JxW[qp];
159 
160  for (unsigned int i=0; i != n_u_dofs; i++)
161  {
162  libMesh::RealGradient u_gradphi( dphi_dxi[i][qp], dphi_deta[i][qp] );
163 
164  for( unsigned int alpha = 0; alpha < dim; alpha++ )
165  {
166  for( unsigned int beta = 0; beta < dim; beta++ )
167  {
168  const libMesh::Real common_factor = _lambda_factor*0.5*this->_h0*jac;
169  const libMesh::Real factor = common_factor*tau(alpha,beta);
170 
171  libMesh::Real u_diagterm = dgradu_dt(beta)*u_gradphi(alpha) + dgradu_dt(alpha)*u_gradphi(beta);
172  libMesh::Real v_diagterm = dgradv_dt(beta)*u_gradphi(alpha) + dgradv_dt(alpha)*u_gradphi(beta);
173  libMesh::Real w_diagterm = dgradw_dt(beta)*u_gradphi(alpha) + dgradw_dt(alpha)*u_gradphi(beta);
174 
175  Fu(i) += factor*u_diagterm;
176  Fv(i) += factor*v_diagterm;
177  Fw(i) += factor*w_diagterm;
178 
179  for( unsigned int lambda = 0; lambda < dim; lambda++ )
180  {
181  for( unsigned int mu = 0; mu < dim; mu++ )
182  {
183  const libMesh::Real C1 = common_factor*C(alpha,beta,lambda,mu);
184 
185  const libMesh::Real gamma_u = 0.5*( dgradu_dt(lambda)*(grad_x(mu)+grad_u(mu)) +
186  (grad_x(lambda)+grad_u(lambda))*dgradu_dt(mu) );
187 
188  const libMesh::Real gamma_v = 0.5*( dgradv_dt(lambda)*(grad_y(mu)+grad_v(mu)) +
189  (grad_y(lambda)+grad_v(lambda))*dgradv_dt(mu) );
190 
191  const libMesh::Real gamma_w = 0.5*( dgradw_dt(lambda)*(grad_z(mu)+grad_w(mu)) +
192  (grad_z(lambda)+grad_w(lambda))*dgradw_dt(mu) );
193 
194  const libMesh::Real x_term = C1*( (grad_x(beta)+grad_u(beta))*u_gradphi(alpha) +
195  (grad_x(alpha)+grad_u(alpha))*u_gradphi(beta) );
196 
197  const libMesh::Real y_term = C1*( (grad_y(beta)+grad_v(beta))*u_gradphi(alpha) +
198  (grad_y(alpha)+grad_v(alpha))*u_gradphi(beta) );
199 
200  const libMesh::Real z_term = C1*( (grad_z(beta)+grad_w(beta))*u_gradphi(alpha) +
201  (grad_z(alpha)+grad_w(alpha))*u_gradphi(beta) );
202 
203  const libMesh::Real gamma_sum = gamma_u + gamma_v + gamma_w;
204 
205  Fu(i) += x_term*gamma_sum;
206 
207  Fv(i) += y_term*gamma_sum;
208 
209  Fw(i) += z_term*gamma_sum;
210  }
211  }
212  }
213  }
214  }
215 
216  if( compute_jacobian )
217  {
218  for (unsigned int i=0; i != n_u_dofs; i++)
219  {
220  libMesh::RealGradient u_gradphi_i( dphi_dxi[i][qp], dphi_deta[i][qp] );
221 
222  for (unsigned int j=0; j != n_u_dofs; j++)
223  {
224  libMesh::RealGradient u_gradphi_j( dphi_dxi[j][qp], dphi_deta[j][qp] );
225 
226  for( unsigned int alpha = 0; alpha < dim; alpha++ )
227  {
228  for( unsigned int beta = 0; beta < dim; beta++ )
229  {
230  const libMesh::Real common_factor = _lambda_factor*0.5*this->_h0*jac;
231 
232  libMesh::Real u_diagterm = dgradu_dt(beta)*u_gradphi_i(alpha) + dgradu_dt(alpha)*u_gradphi_i(beta);
233  libMesh::Real v_diagterm = dgradv_dt(beta)*u_gradphi_i(alpha) + dgradv_dt(alpha)*u_gradphi_i(beta);
234  libMesh::Real w_diagterm = dgradw_dt(beta)*u_gradphi_i(alpha) + dgradw_dt(alpha)*u_gradphi_i(beta);
235 
236  const libMesh::Real factor = common_factor*tau(alpha,beta);
237 
238  const libMesh::Real ddiagterm = u_gradphi_j(beta)*u_gradphi_i(alpha) + u_gradphi_j(alpha)*u_gradphi_i(beta);
239 
240  Kuu(i,j) += factor*ddiagterm*context.get_elem_solution_rate_derivative();
241  Kvv(i,j) += factor*ddiagterm*context.get_elem_solution_rate_derivative();
242  Kww(i,j) += factor*ddiagterm*context.get_elem_solution_rate_derivative();
243 
244  for( unsigned int lambda = 0; lambda < dim; lambda++ )
245  {
246  for( unsigned int mu = 0; mu < dim; mu++ )
247  {
248  const libMesh::Real dgamma_du = 0.5*( u_gradphi_j(lambda)*(grad_x(mu)+grad_u(mu)) +
249  (grad_x(lambda)+grad_u(lambda))*u_gradphi_j(mu) );
250 
251  const libMesh::Real dgamma_dv = 0.5*( u_gradphi_j(lambda)*(grad_y(mu)+grad_v(mu)) +
252  (grad_y(lambda)+grad_v(lambda))*u_gradphi_j(mu) );
253 
254  const libMesh::Real dgamma_dw = 0.5*( u_gradphi_j(lambda)*(grad_z(mu)+grad_w(mu)) +
255  (grad_z(lambda)+grad_w(lambda))*u_gradphi_j(mu) );
256 
257  const libMesh::Real C1 = common_factor*C(alpha,beta,lambda,mu);
258 
259  const libMesh::Real dfactor_du = C1*dgamma_du*context.get_elem_solution_derivative();
260  const libMesh::Real dfactor_dv = C1*dgamma_dv*context.get_elem_solution_derivative();
261  const libMesh::Real dfactor_dw = C1*dgamma_dw*context.get_elem_solution_derivative();
262 
263  Kuu(i,j) += u_diagterm*dfactor_du;
264  Kuv(i,j) += u_diagterm*dfactor_dv;
265  Kuw(i,j) += u_diagterm*dfactor_dw;
266 
267  Kvu(i,j) += v_diagterm*dfactor_du;
268  Kvv(i,j) += v_diagterm*dfactor_dv;
269  Kvw(i,j) += v_diagterm*dfactor_dw;
270 
271  Kwu(i,j) += w_diagterm*dfactor_du;
272  Kwv(i,j) += w_diagterm*dfactor_dv;
273  Kww(i,j) += w_diagterm*dfactor_dw;
274 
275  const libMesh::Real gamma_u = 0.5*( dgradu_dt(lambda)*(grad_x(mu)+grad_u(mu)) +
276  (grad_x(lambda)+grad_u(lambda))*dgradu_dt(mu) );
277 
278  const libMesh::Real gamma_v = 0.5*( dgradv_dt(lambda)*(grad_y(mu)+grad_v(mu)) +
279  (grad_y(lambda)+grad_v(lambda))*dgradv_dt(mu) );
280 
281  const libMesh::Real gamma_w = 0.5*( dgradw_dt(lambda)*(grad_z(mu)+grad_w(mu)) +
282  (grad_z(lambda)+grad_w(lambda))*dgradw_dt(mu) );
283 
284  const libMesh::Real x_term = C1*( (grad_x(beta)+grad_u(beta))*u_gradphi_i(alpha) +
285  (grad_x(alpha)+grad_u(alpha))*u_gradphi_i(beta) );
286 
287  const libMesh::Real y_term = C1*( (grad_y(beta)+grad_v(beta))*u_gradphi_i(alpha) +
288  (grad_y(alpha)+grad_v(alpha))*u_gradphi_i(beta) );
289 
290  const libMesh::Real z_term = C1*( (grad_z(beta)+grad_w(beta))*u_gradphi_i(alpha) +
291  (grad_z(alpha)+grad_w(alpha))*u_gradphi_i(beta) );
292 
293  const libMesh::Real gamma_sum = gamma_u + gamma_v + gamma_w;
294 
295  // Here, we're missing derivatives of C(alpha,beta,lambda,mu) w.r.t. strain
296  // Nonzero for hyperelasticity models
297  const libMesh::Real dxterm_du = C1*( u_gradphi_j(beta)*u_gradphi_i(alpha)
298  + u_gradphi_j(alpha)*u_gradphi_i(beta) )*context.get_elem_solution_derivative();
299 
300  const libMesh::Real dyterm_dv = dxterm_du;
301  const libMesh::Real dzterm_dw = dxterm_du;
302 
303  Kuu(i,j) += gamma_sum*dxterm_du;
304 
305  Kvv(i,j) += gamma_sum*dyterm_dv;
306 
307  Kww(i,j) += gamma_sum*dzterm_dw;
308 
309  const libMesh::Real dgamma_sum_du = 0.5*( u_gradphi_j(lambda)*(grad_x(mu)+grad_u(mu))*context.get_elem_solution_rate_derivative()
310  + dgradu_dt(lambda)*u_gradphi_j(mu)*context.get_elem_solution_derivative()
311  + u_gradphi_j(mu)*(grad_x(lambda)+grad_u(lambda))*context.get_elem_solution_rate_derivative()
312  + dgradu_dt(mu)*u_gradphi_j(lambda)*context.get_elem_solution_derivative() );
313 
314  const libMesh::Real dgamma_sum_dv = 0.5*( u_gradphi_j(lambda)*(grad_y(mu)+grad_v(mu))*context.get_elem_solution_rate_derivative()
315  + dgradv_dt(lambda)*u_gradphi_j(mu)*context.get_elem_solution_derivative()
316  + u_gradphi_j(mu)*(grad_y(lambda)+grad_v(lambda))*context.get_elem_solution_rate_derivative()
317  + dgradv_dt(mu)*u_gradphi_j(lambda)*context.get_elem_solution_derivative() );
318 
319  const libMesh::Real dgamma_sum_dw = 0.5*( u_gradphi_j(lambda)*(grad_z(mu)+grad_w(mu))*context.get_elem_solution_rate_derivative()
320  + dgradw_dt(lambda)*u_gradphi_j(mu)*context.get_elem_solution_derivative()
321  + u_gradphi_j(mu)*(grad_z(lambda)+grad_w(lambda))*context.get_elem_solution_rate_derivative()
322  + dgradw_dt(mu)*u_gradphi_j(lambda)*context.get_elem_solution_derivative() );
323 
324  Kuu(i,j) += x_term*dgamma_sum_du;
325  Kuv(i,j) += x_term*dgamma_sum_dv;
326  Kuw(i,j) += x_term*dgamma_sum_dw;
327 
328  Kvu(i,j) += y_term*dgamma_sum_du;
329  Kvv(i,j) += y_term*dgamma_sum_dv;
330  Kvw(i,j) += y_term*dgamma_sum_dw;
331 
332  Kwu(i,j) += z_term*dgamma_sum_du;
333  Kwv(i,j) += z_term*dgamma_sum_dv;
334  Kww(i,j) += z_term*dgamma_sum_dw;
335  }
336  }
337  }
338  }
339  }
340  }
341  }
342 
343  }
344  }
345 
346 
347 } // end namespace GRINS
GRINS namespace.
static PhysicsName elastic_membrane()
std::string PhysicsName
virtual void damping_residual(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
void parse_enabled_subdomains(const GetPot &input, const std::string &physics_name)
Definition: physics.C:64
static PhysicsName elastic_membrane_rayleigh_damping()

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