GRINS-0.7.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-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_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 
66  AssemblyContext& context,
67  CachedValues& /*cache*/ )
68  {
69  const unsigned int n_u_dofs = context.get_dof_indices(_disp_vars.u()).size();
70 
71  const std::vector<libMesh::Real> &JxW =
72  this->get_fe(context)->get_JxW();
73 
74  const std::vector<std::vector<libMesh::Real> >& u_phi =
75  this->get_fe(context)->get_phi();
76 
77  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(_disp_vars.u());
78  libMesh::DenseSubVector<libMesh::Number> &Fv = context.get_elem_residual(_disp_vars.v());
79  libMesh::DenseSubVector<libMesh::Number> &Fw = context.get_elem_residual(_disp_vars.w());
80 
81  libMesh::DenseSubMatrix<libMesh::Number>& Kuv = context.get_elem_jacobian(_disp_vars.u(),_disp_vars.v());
82  libMesh::DenseSubMatrix<libMesh::Number>& Kuw = context.get_elem_jacobian(_disp_vars.u(),_disp_vars.w());
83 
84  libMesh::DenseSubMatrix<libMesh::Number>& Kvu = context.get_elem_jacobian(_disp_vars.v(),_disp_vars.u());
85  libMesh::DenseSubMatrix<libMesh::Number>& Kvw = context.get_elem_jacobian(_disp_vars.v(),_disp_vars.w());
86 
87  libMesh::DenseSubMatrix<libMesh::Number>& Kwu = context.get_elem_jacobian(_disp_vars.w(),_disp_vars.u());
88  libMesh::DenseSubMatrix<libMesh::Number>& Kwv = context.get_elem_jacobian(_disp_vars.w(),_disp_vars.v());
89 
90  unsigned int n_qpoints = context.get_element_qrule().n_points();
91 
92  // All shape function gradients are w.r.t. master element coordinates
93  const std::vector<std::vector<libMesh::Real> >& dphi_dxi =
94  this->get_fe(context)->get_dphidxi();
95 
96  const std::vector<std::vector<libMesh::Real> >& dphi_deta =
97  this->get_fe(context)->get_dphideta();
98 
99  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( _disp_vars.u() );
100  const libMesh::DenseSubVector<libMesh::Number>& v_coeffs = context.get_elem_solution( _disp_vars.v() );
101  const libMesh::DenseSubVector<libMesh::Number>& w_coeffs = context.get_elem_solution( _disp_vars.w() );
102 
103  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
104  const std::vector<libMesh::RealGradient>& dxdeta = this->get_fe(context)->get_dxyzdeta();
105 
106  for (unsigned int qp=0; qp != n_qpoints; qp++)
107  {
108  // sqrt(det(a_cov)), a_cov being the covariant metric tensor of undeformed body
109  libMesh::Real sqrt_a = sqrt( dxdxi[qp]*dxdxi[qp]*dxdeta[qp]*dxdeta[qp]
110  - dxdxi[qp]*dxdeta[qp]*dxdeta[qp]*dxdxi[qp] );
111 
112  // Gradients are w.r.t. master element coordinates
113  libMesh::Gradient grad_u, grad_v, grad_w;
114  for( unsigned int d = 0; d < n_u_dofs; d++ )
115  {
116  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp], dphi_deta[d][qp] );
117  grad_u += u_coeffs(d)*u_gradphi;
118  grad_v += v_coeffs(d)*u_gradphi;
119  grad_w += w_coeffs(d)*u_gradphi;
120  }
121 
122  libMesh::RealGradient dudxi( grad_u(0), grad_v(0), grad_w(0) );
123  libMesh::RealGradient dudeta( grad_u(1), grad_v(1), grad_w(1) );
124 
125  libMesh::RealGradient A_1 = dxdxi[qp] + dudxi;
126  libMesh::RealGradient A_2 = dxdeta[qp] + dudeta;
127 
128  libMesh::RealGradient A_3 = A_1.cross(A_2);
129 
130  /* The formula here is actually
131  P*\sqrt{\frac{A}{a}}*A_3, where A_3 is a unit vector
132  But, |A_3| = \sqrt{A} so the normalizing part kills
133  the \sqrt{A} in the numerator, so we can leave it out
134  and *not* normalize A_3.
135  */
136  libMesh::RealGradient traction = _pressure/sqrt_a*A_3;
137 
138  libMesh::Real jac = JxW[qp];
139 
140  for (unsigned int i=0; i != n_u_dofs; i++)
141  {
142  Fu(i) -= traction(0)*u_phi[i][qp]*jac;
143 
144  Fv(i) -= traction(1)*u_phi[i][qp]*jac;
145 
146  Fw(i) -= traction(2)*u_phi[i][qp]*jac;
147 
148  if( compute_jacobian )
149  {
150  for (unsigned int j=0; j != n_u_dofs; j++)
151  {
152  libMesh::RealGradient u_gradphi( dphi_dxi[j][qp], dphi_deta[j][qp] );
153 
154  const libMesh::Real dt0_dv = _pressure/sqrt_a*(u_gradphi(0)*A_2(2) - A_1(2)*u_gradphi(1));
155  const libMesh::Real dt0_dw = _pressure/sqrt_a*(A_1(1)*u_gradphi(1) - u_gradphi(0)*A_2(1));
156 
157  const libMesh::Real dt1_du = _pressure/sqrt_a*(A_1(2)*u_gradphi(1) - u_gradphi(0)*A_2(2));
158  const libMesh::Real dt1_dw = _pressure/sqrt_a*(u_gradphi(0)*A_2(0) - A_1(0)*u_gradphi(1));
159 
160  const libMesh::Real dt2_du = _pressure/sqrt_a*(u_gradphi(0)*A_2(1) - A_1(1)*u_gradphi(1));
161  const libMesh::Real dt2_dv = _pressure/sqrt_a*(A_1(0)*u_gradphi(1) - u_gradphi(0)*A_2(0));
162 
163  Kuv(i,j) -= dt0_dv*u_phi[i][qp]*jac;
164  Kuw(i,j) -= dt0_dw*u_phi[i][qp]*jac;
165 
166  Kvu(i,j) -= dt1_du*u_phi[i][qp]*jac;
167  Kvw(i,j) -= dt1_dw*u_phi[i][qp]*jac;
168 
169  Kwu(i,j) -= dt2_du*u_phi[i][qp]*jac;
170  Kwv(i,j) -= dt2_dv*u_phi[i][qp]*jac;
171  }
172  }
173  }
174  }
175 
176  return;
177  }
178 
179  void ElasticMembraneConstantPressure::reset_pressure( libMesh::Real pressure_in )
180  {
181  _pressure = pressure_in;
182  return;
183  }
184 
185 } // 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.
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context, CachedValues &cache)
Time dependent part(s) of physics for element interiors.
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:64
const libMesh::FEGenericBase< libMesh::Real > * get_fe(const AssemblyContext &context)

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