GRINS-0.8.0
elastic_membrane_constant_pressure.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // GRINS - General Reacting Incompressible Navier-Stokes
5 //
6 // Copyright (C) 2014-2017 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_config.h"
30 #include "grins/assembly_context.h"
31 
32 // libMesh
33 #include "libmesh/getpot.h"
34 #include "libmesh/quadrature.h"
35 #include "libmesh/boundary_info.h"
36 #include "libmesh/fem_system.h"
37 
38 namespace GRINS
39 {
41  : ElasticMembraneAbstract(physics_name,input),
42  _pressure(0.0)
43  {
44  if( !input.have_variable("Physics/"+physics_name+"/pressure") )
45  {
46  std::cerr << "Error: Must input pressure for "+physics_name+"." << std::endl
47  << " Please set Physics/"+physics_name+"/pressure." << std::endl;
48  libmesh_error();
49  }
50 
51  this->set_parameter
52  (_pressure, input, "Physics/"+physics_name+"/pressure", _pressure );
53 
54  // If the user specified enabled subdomains in this Physics section,
55  // that's an error; we're slave to ElasticMembrane.
56  if( input.have_variable("Physics/"+PhysicsNaming::elastic_membrane_constant_pressure()+"/enabled_subdomains" ) )
57  libmesh_error_msg("ERROR: Cannot specify subdomains for "
59  +"! Must specify subdomains through "
61 
63 
64  if( this->_disp_vars.dim() != 3 )
65  libmesh_error_msg("ERROR: ElasticMembraneConstantPressure only valid for three dimensions! Make sure you have three components in your Displacement type variable.");
66  }
67 
69  ( bool compute_jacobian,
70  AssemblyContext & context )
71  {
72  const unsigned int n_u_dofs = context.get_dof_indices(_disp_vars.u()).size();
73 
74  const std::vector<libMesh::Real> &JxW =
75  this->get_fe(context)->get_JxW();
76 
77  const std::vector<std::vector<libMesh::Real> >& u_phi =
78  this->get_fe(context)->get_phi();
79 
80  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_disp_vars.u());
81  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_disp_vars.v());
82  libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(_disp_vars.w());
83 
84  libMesh::DenseSubMatrix<libMesh::Number>& Kuv = context.get_elem_jacobian(_disp_vars.u(),_disp_vars.v());
85  libMesh::DenseSubMatrix<libMesh::Number>& Kuw = context.get_elem_jacobian(_disp_vars.u(),_disp_vars.w());
86 
87  libMesh::DenseSubMatrix<libMesh::Number>& Kvu = context.get_elem_jacobian(_disp_vars.v(),_disp_vars.u());
88  libMesh::DenseSubMatrix<libMesh::Number>& Kvw = context.get_elem_jacobian(_disp_vars.v(),_disp_vars.w());
89 
90  libMesh::DenseSubMatrix<libMesh::Number>& Kwu = context.get_elem_jacobian(_disp_vars.w(),_disp_vars.u());
91  libMesh::DenseSubMatrix<libMesh::Number>& Kwv = context.get_elem_jacobian(_disp_vars.w(),_disp_vars.v());
92 
93  unsigned int n_qpoints = context.get_element_qrule().n_points();
94 
95  // All shape function gradients are w.r.t. master element coordinates
96  const std::vector<std::vector<libMesh::Real> >& dphi_dxi =
97  this->get_fe(context)->get_dphidxi();
98 
99  const std::vector<std::vector<libMesh::Real> >& dphi_deta =
100  this->get_fe(context)->get_dphideta();
101 
102  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( _disp_vars.u() );
103  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( _disp_vars.v() );
104  const libMesh::DenseSubVector<libMesh::Number>& w_coeffs = context.get_elem_solution( _disp_vars.w() );
105 
106  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
107  const std::vector<libMesh::RealGradient>& dxdeta = this->get_fe(context)->get_dxyzdeta();
108 
109  for (unsigned int qp=0; qp != n_qpoints; qp++)
110  {
111  // sqrt(det(a_cov)), a_cov being the covariant metric tensor of undeformed body
112  libMesh::Real sqrt_a = sqrt( dxdxi[qp]*dxdxi[qp]*dxdeta[qp]*dxdeta[qp]
113  - dxdxi[qp]*dxdeta[qp]*dxdeta[qp]*dxdxi[qp] );
114 
115  // Gradients are w.r.t. master element coordinates
116  libMesh::Gradient grad_u, grad_v, grad_w;
117  for( unsigned int d = 0; d < n_u_dofs; d++ )
118  {
119  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp], dphi_deta[d][qp] );
120  grad_u += u_coeffs(d)*u_gradphi;
121  grad_v += v_coeffs(d)*u_gradphi;
122  grad_w += w_coeffs(d)*u_gradphi;
123  }
124 
125  libMesh::RealGradient dudxi( grad_u(0), grad_v(0), grad_w(0) );
126  libMesh::RealGradient dudeta( grad_u(1), grad_v(1), grad_w(1) );
127 
128  libMesh::RealGradient A_1 = dxdxi[qp] + dudxi;
129  libMesh::RealGradient A_2 = dxdeta[qp] + dudeta;
130 
131  libMesh::RealGradient A_3 = A_1.cross(A_2);
132 
133  /* The formula here is actually
134  P*\sqrt{\frac{A}{a}}*A_3, where A_3 is a unit vector
135  But, |A_3| = \sqrt{A} so the normalizing part kills
136  the \sqrt{A} in the numerator, so we can leave it out
137  and *not* normalize A_3.
138  */
139  libMesh::RealGradient traction = _pressure/sqrt_a*A_3;
140 
141  libMesh::Real jac = JxW[qp];
142 
143  for (unsigned int i=0; i != n_u_dofs; i++)
144  {
145  Fu(i) -= traction(0)*u_phi[i][qp]*jac;
146 
147  Fv(i) -= traction(1)*u_phi[i][qp]*jac;
148 
149  Fw(i) -= traction(2)*u_phi[i][qp]*jac;
150 
151  if( compute_jacobian )
152  {
153  for (unsigned int j=0; j != n_u_dofs; j++)
154  {
155  libMesh::RealGradient u_gradphi( dphi_dxi[j][qp], dphi_deta[j][qp] );
156 
157  const libMesh::Real dt0_dv = _pressure/sqrt_a*(u_gradphi(0)*A_2(2) - A_1(2)*u_gradphi(1));
158  const libMesh::Real dt0_dw = _pressure/sqrt_a*(A_1(1)*u_gradphi(1) - u_gradphi(0)*A_2(1));
159 
160  const libMesh::Real dt1_du = _pressure/sqrt_a*(A_1(2)*u_gradphi(1) - u_gradphi(0)*A_2(2));
161  const libMesh::Real dt1_dw = _pressure/sqrt_a*(u_gradphi(0)*A_2(0) - A_1(0)*u_gradphi(1));
162 
163  const libMesh::Real dt2_du = _pressure/sqrt_a*(u_gradphi(0)*A_2(1) - A_1(1)*u_gradphi(1));
164  const libMesh::Real dt2_dv = _pressure/sqrt_a*(A_1(0)*u_gradphi(1) - u_gradphi(0)*A_2(0));
165 
166  Kuv(i,j) -= dt0_dv*u_phi[i][qp]*jac;
167  Kuw(i,j) -= dt0_dw*u_phi[i][qp]*jac;
168 
169  Kvu(i,j) -= dt1_du*u_phi[i][qp]*jac;
170  Kvw(i,j) -= dt1_dw*u_phi[i][qp]*jac;
171 
172  Kwu(i,j) -= dt2_du*u_phi[i][qp]*jac;
173  Kwv(i,j) -= dt2_dv*u_phi[i][qp]*jac;
174  }
175  }
176  }
177  }
178 
179  return;
180  }
181 
182  void ElasticMembraneConstantPressure::reset_pressure( libMesh::Real pressure_in )
183  {
184  _pressure = pressure_in;
185  return;
186  }
187 
188 } // end namespace GRINS
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.
static PhysicsName elastic_membrane_constant_pressure()
GRINS namespace.
static PhysicsName elastic_membrane()
std::string PhysicsName
void parse_enabled_subdomains(const GetPot &input, const std::string &physics_name)
Definition: physics.C:65
unsigned int dim() const
Number of components.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.

Generated on Tue Dec 19 2017 12:47:28 for GRINS-0.8.0 by  doxygen 1.8.9.1