GRINS-0.8.0
elastic_cable.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
26 #include "grins/elastic_cable.h"
27 
28 // GRINS
29 #include "grins_config.h"
33 
34 // libMesh
35 #include "libmesh/getpot.h"
36 #include "libmesh/quadrature.h"
37 #include "libmesh/boundary_info.h"
38 #include "libmesh/fem_system.h"
39 #include "libmesh/elem.h"
40 
41 namespace GRINS
42 {
43  template<typename StressStrainLaw>
44  ElasticCable<StressStrainLaw>::ElasticCable( const PhysicsName& physics_name, const GetPot& input,
45  bool is_compressible )
46  : ElasticCableBase<StressStrainLaw>(physics_name,input,is_compressible)
47  {
48  this->_ic_handler = new GenericICHandler(physics_name, input);
49  }
50 
51  template<typename StressStrainLaw>
54  {
55  std::string section = "Physics/"+PhysicsNaming::elastic_cable()+"/output_vars";
56 
57  if( input.have_variable(section) )
58  {
59  unsigned int n_vars = input.vector_variable_size(section);
60 
61  for( unsigned int v = 0; v < n_vars; v++ )
62  {
63  std::string name = input(section,"DIE!",v);
64 
65  if( name == std::string("stress") )
66  {
67  // sigma_xx only, i.e., locally normal to the cutting plane
68  // sigma_yy=sigma_zz = 0 by assumption of this Physics
69  _stress_indices.resize(1);
70 
71  this->_stress_indices[0] = postprocessing.register_quantity("cable_stress");
72 
73  }
74  else if( name == std::string("strain") )
75  {
76  // eps_xx only, i.e., locally normal to the cutting plane
77  // eps_yy=eps_zz = 0 by assumption of this Physics
78  _strain_indices.resize(1);
79 
80  this->_strain_indices[0] = postprocessing.register_quantity("cable_strain");
81 
82  }
83  else if( name == std::string("force") )
84  {
85  // force_x only, i.e., locally normal to the cutting plane
86  // force_y=force_z=0 by assumption of this Physics
87  _force_indices.resize(1);
88 
89  this->_force_indices[0] = postprocessing.register_quantity("cable_force");
90 
91  }
92  else
93  {
94  std::cerr << "Error: Invalue output_vars value for "+PhysicsNaming::elastic_cable() << std::endl
95  << " Found " << name << std::endl
96  << " Acceptable values are: stress" << std::endl
97  << " strain" << std::endl
98  << " force " << std::endl;
99  libmesh_error();
100  }
101  }
102  }
103  return;
104  }
105 
106  template<typename StressStrainLaw>
108  const AssemblyContext& context,
109  const libMesh::Point& point,
110  libMesh::Real& value )
111  {
112  bool is_strain = false;
113  if( !_strain_indices.empty() )
114  is_strain = ( _strain_indices[0] == quantity_index );
115 
116  bool is_stress = false;
117  if( !_stress_indices.empty() )
118  is_stress = ( _stress_indices[0] == quantity_index );
119 
120  bool is_force = false;
121  if( !_force_indices.empty() )
122  is_force = ( _force_indices[0] == quantity_index );
123 
124  if( is_strain || is_stress || is_force )
125  {
126  const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
127 
128  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
129  const libMesh::DenseSubVector<libMesh::Number>* v_coeffs = NULL;
130 
131  const libMesh::DenseSubVector<libMesh::Number>* w_coeffs = NULL;
132 
133  if( this->_disp_vars.dim() >= 2 )
134  v_coeffs = &context.get_elem_solution( this->_disp_vars.v() );
135 
136  if( this->_disp_vars.dim() == 3 )
137  w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
138 
139  // Build new FE for the current point. We need this to build tensors at point.
140  libMesh::UniquePtr<libMesh::FEGenericBase<libMesh::Real> > fe_new = this->build_new_fe( &context.get_elem(), this->get_fe(context), point );
141 
142  const std::vector<std::vector<libMesh::Real> >& dphi_dxi = fe_new->get_dphidxi();
143 
144  // Need these to build up the covariant and contravariant metric tensors
145  const std::vector<libMesh::RealGradient>& dxdxi = fe_new->get_dxyzdxi();
146 
147  // Gradients are w.r.t. master element coordinates
148  libMesh::Gradient grad_u, grad_v, grad_w;
149  for( unsigned int d = 0; d < n_u_dofs; d++ )
150  {
151  libMesh::RealGradient u_gradphi( dphi_dxi[d][0] );
152  grad_u += u_coeffs(d)*u_gradphi;
153 
154  if( this->_disp_vars.dim() >= 2 )
155  grad_v += (*v_coeffs)(d)*u_gradphi;
156 
157  if( this->_disp_vars.dim() == 3 )
158  grad_w += (*w_coeffs)(d)*u_gradphi;
159  }
160 
161  libMesh::RealGradient grad_x( dxdxi[0](0) );
162  libMesh::RealGradient grad_y( dxdxi[0](1) );
163  libMesh::RealGradient grad_z( dxdxi[0](2) );
164 
165  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
166  libMesh::Real lambda_sq = 0;
167 
168  this->compute_metric_tensors(0, *fe_new, context,
169  grad_u, grad_v, grad_w,
170  a_cov, a_contra, A_cov, A_contra, lambda_sq );
171 
172  libMesh::Real det_a = a_cov(0,0)*a_cov(1,1) - a_cov(0,1)*a_cov(1,0);
173  libMesh::Real det_A = A_cov(0,0)*A_cov(1,1) - A_cov(0,1)*A_cov(1,0);
174 
175  libMesh::Real I3 = lambda_sq*det_A/det_a;
176 
178  this->_stress_strain_law.compute_stress(1,a_contra,a_cov,A_contra,A_cov,tau);
179 
180  // We have everything we need for strain now, so check if we are computing strain
181  if( is_strain )
182  {
183  if( _strain_indices[0] == quantity_index )
184  {
185  value = 0.5*(A_cov(0,0) - a_cov(0,0));
186  }
187  else
188  {
189  //Wat?!
190  libmesh_error();
191  }
192  return;
193  }
194 
195  if( is_stress )
196  {
197  if( _stress_indices[0] == quantity_index )
198  {
199  // Need to convert to Cauchy stress
200  value = tau(0,0)/std::sqrt(I3);
201  }
202  else
203  {
204  libmesh_error();
205  }
206  }
207  if( is_force )
208  {
209 
210  if( _force_indices[0] == quantity_index )
211  {
212  //This force is in deformed configuration and will be significantly influenced by large strains
213  value = tau(0,0)/std::sqrt(I3)*this->_A;
214  }
215  else
216  {
217  libmesh_error();
218  }
219  }
220  }
221  }
222 
223  template<typename StressStrainLaw>
225  ( bool compute_jacobian, AssemblyContext & context )
226  {
227  const unsigned int n_u_dofs = context.get_dof_indices(this->_disp_vars.u()).size();
228 
229  const std::vector<libMesh::Real> &JxW =
230  this->get_fe(context)->get_JxW();
231 
232  // Residuals that we're populating
233  libMesh::DenseSubVector<libMesh::Number> &Fu = context.get_elem_residual(this->_disp_vars.u());
234  libMesh::DenseSubVector<libMesh::Number>* Fv = NULL;
235  libMesh::DenseSubVector<libMesh::Number>* Fw = NULL;
236 
237  libMesh::DenseSubMatrix<libMesh::Number>& Kuu = context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.u());
238  libMesh::DenseSubMatrix<libMesh::Number>* Kuv = NULL;
239  libMesh::DenseSubMatrix<libMesh::Number>* Kuw = NULL;
240 
241  libMesh::DenseSubMatrix<libMesh::Number>* Kvu = NULL;
242  libMesh::DenseSubMatrix<libMesh::Number>* Kvv = NULL;
243  libMesh::DenseSubMatrix<libMesh::Number>* Kvw = NULL;
244 
245  libMesh::DenseSubMatrix<libMesh::Number>* Kwu = NULL;
246  libMesh::DenseSubMatrix<libMesh::Number>* Kwv = NULL;
247  libMesh::DenseSubMatrix<libMesh::Number>* Kww = NULL;
248 
249  if( this->_disp_vars.dim() >= 2 )
250  {
251  Fv = &context.get_elem_residual(this->_disp_vars.v());
252 
253  Kuv = &context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.v());
254  Kvu = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.u());
255  Kvv = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.v());
256  }
257 
258  if( this->_disp_vars.dim() == 3 )
259  {
260  Fw = &context.get_elem_residual(this->_disp_vars.w());
261 
262  Kuw = &context.get_elem_jacobian(this->_disp_vars.u(),this->_disp_vars.w());
263  Kvw = &context.get_elem_jacobian(this->_disp_vars.v(),this->_disp_vars.w());
264  Kwu = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.u());
265  Kwv = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.v());
266  Kww = &context.get_elem_jacobian(this->_disp_vars.w(),this->_disp_vars.w());
267  }
268 
269 
270  unsigned int n_qpoints = context.get_element_qrule().n_points();
271 
272  // All shape function gradients are w.r.t. master element coordinates
273  const std::vector<std::vector<libMesh::Real> >& dphi_dxi = this->get_fe(context)->get_dphidxi();
274 
275  const libMesh::DenseSubVector<libMesh::Number>& u_coeffs = context.get_elem_solution( this->_disp_vars.u() );
276  const libMesh::DenseSubVector<libMesh::Number>* v_coeffs = NULL;
277  const libMesh::DenseSubVector<libMesh::Number>* w_coeffs = NULL;
278 
279  if( this->_disp_vars.dim() >= 2 )
280  v_coeffs = &context.get_elem_solution( this->_disp_vars.v() );
281 
282  if( this->_disp_vars.dim() == 3 )
283  w_coeffs = &context.get_elem_solution( this->_disp_vars.w() );
284 
285  // Need these to build up the covariant and contravariant metric tensors
286  const std::vector<libMesh::RealGradient>& dxdxi = this->get_fe(context)->get_dxyzdxi();
287 
288  const unsigned int dim = 1; // The cable dimension is always 1 for this physics
289 
290  for (unsigned int qp=0; qp != n_qpoints; qp++)
291  {
292  // Gradients are w.r.t. master element coordinates
293  libMesh::Gradient grad_u, grad_v, grad_w;
294 
295  for( unsigned int d = 0; d < n_u_dofs; d++ )
296  {
297  libMesh::RealGradient u_gradphi( dphi_dxi[d][qp] );
298  grad_u += u_coeffs(d)*u_gradphi;
299 
300  if( this->_disp_vars.dim() >= 2 )
301  grad_v += (*v_coeffs)(d)*u_gradphi;
302 
303  if( this->_disp_vars.dim() == 3 )
304  grad_w += (*w_coeffs)(d)*u_gradphi;
305  }
306 
307  libMesh::RealGradient grad_x( dxdxi[qp](0) );
308  libMesh::RealGradient grad_y( dxdxi[qp](1) );
309  libMesh::RealGradient grad_z( dxdxi[qp](2) );
310 
311  libMesh::TensorValue<libMesh::Real> a_cov, a_contra, A_cov, A_contra;
312  libMesh::Real lambda_sq = 0;
313 
314  this->compute_metric_tensors( qp, *(this->get_fe(context)), context,
315  grad_u, grad_v, grad_w,
316  a_cov, a_contra, A_cov, A_contra,
317  lambda_sq );
318 
319  // Compute stress tensor
322  this->_stress_strain_law.compute_stress_and_elasticity(dim,a_contra,a_cov,A_contra,A_cov,tau,C);
323 
324 
325  libMesh::Real jac = JxW[qp];
326 
327  for (unsigned int i=0; i != n_u_dofs; i++)
328  {
329  libMesh::RealGradient u_gradphi( dphi_dxi[i][qp] );
330 
331  const libMesh::Real res_term = tau(0,0)*this->_A*jac*u_gradphi(0);
332 
333  Fu(i) += res_term*(grad_x(0) + grad_u(0));
334 
335  if( this->_disp_vars.dim() >= 2 )
336  (*Fv)(i) += res_term*(grad_y(0) + grad_v(0));
337 
338  if( this->_disp_vars.dim() == 3 )
339  (*Fw)(i) += res_term*(grad_z(0) + grad_w(0));
340  }
341 
342  if( compute_jacobian )
343  {
344  for(unsigned int i=0; i != n_u_dofs; i++)
345  {
346  libMesh::RealGradient u_gradphi_I( dphi_dxi[i][qp] );
347  for(unsigned int j=0; j != n_u_dofs; j++)
348  {
349  libMesh::RealGradient u_gradphi_J( dphi_dxi[j][qp] );
350 
351  const libMesh::Real diag_term = this->_A*jac*tau(0,0)*( u_gradphi_J(0)*u_gradphi_I(0))*context.get_elem_solution_derivative();
352 
353  Kuu(i,j) += diag_term;
354 
355  if( this->_disp_vars.dim() >= 2 )
356  (*Kvv)(i,j) += diag_term;
357 
358  if( this->_disp_vars.dim() == 3 )
359  (*Kww)(i,j) += diag_term;
360 
361  const libMesh::Real dgamma_du = ( u_gradphi_J(0)*(grad_x(0)+grad_u(0)) );
362 
363  const libMesh::Real C1 = this->_A*jac*C(0,0,0,0)*context.get_elem_solution_derivative();
364 
365  const libMesh::Real x_term = C1*( (grad_x(0)+grad_u(0))*u_gradphi_I(0) );
366 
367  Kuu(i,j) += x_term*dgamma_du;
368 
369  libMesh::Real y_term = 0.0;
370  libMesh::Real dgamma_dv = 0.0;
371 
372  if( this->_disp_vars.dim() >= 2 )
373  {
374  dgamma_dv = ( u_gradphi_J(0)*(grad_y(0)+grad_v(0)) );
375  y_term = C1*( (grad_y(0)+grad_v(0))*u_gradphi_I(0) );
376 
377  (*Kuv)(i,j) += x_term*dgamma_dv;
378  (*Kvu)(i,j) += y_term*dgamma_du;
379  (*Kvv)(i,j) += y_term*dgamma_dv;
380  }
381 
382  libMesh::Real z_term = 0.0;
383 
384  if( this->_disp_vars.dim() == 3 )
385  {
386  const libMesh::Real dgamma_dw = ( u_gradphi_J(0)*(grad_z(0)+grad_w(0)) );
387  z_term = C1*( (grad_z(0)+grad_w(0))*u_gradphi_I(0) );
388 
389  (*Kuw)(i,j) += x_term*dgamma_dw;
390  (*Kvw)(i,j) += y_term*dgamma_dw;
391  (*Kwu)(i,j) += z_term*dgamma_du;
392  (*Kwv)(i,j) += z_term*dgamma_dv;
393  (*Kww)(i,j) += z_term*dgamma_dw;
394  }
395 
396  } // end j-loop
397  } // end i-loop
398  } // end if(compute_jacobian)
399  } // end qp loop
400  }
401 
402 } // end namespace GRINS
unsigned int register_quantity(std::string name)
Register quantity to be postprocessed.
GRINS::ICHandlingBase * _ic_handler
Definition: physics.h:258
virtual void compute_postprocessed_quantity(unsigned int quantity_index, const AssemblyContext &context, const libMesh::Point &point, libMesh::Real &value)
Compute the registered postprocessed quantities.
static PhysicsName elastic_cable()
Base class for reading and handling initial conditions for physics classes.
GRINS namespace.
virtual void register_postprocessing_vars(const GetPot &input, PostProcessedQuantities< libMesh::Real > &postprocessing)
Register postprocessing variables for ElasticCable.
Definition: elastic_cable.C:52
virtual void element_time_derivative(bool compute_jacobian, AssemblyContext &context)
Time dependent part(s) of physics for element interiors.
std::string PhysicsName

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