time_harmonic_linear_elasticity_elements.cc
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2024 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Non-inline functions for elements that solve the equations of linear
27 // elasticity in cartesian coordinates
28 
30 
31 
32 namespace oomph
33 {
34  /// Static default value for square of frequency
35  template<unsigned DIM>
36  double
38  1.0;
39 
40 
41  /// ///////////////////////////////////////////////////////////////
42  /// ///////////////////////////////////////////////////////////////
43  /// ///////////////////////////////////////////////////////////////
44 
45  //======================================================================
46  /// Compute the strain tensor at local coordinate s
47  //======================================================================
48  template<unsigned DIM>
50  const Vector<double>& s, DenseMatrix<std::complex<double>>& strain) const
51  {
52 #ifdef PARANOID
53  if ((strain.ncol() != DIM) || (strain.nrow() != DIM))
54  {
55  std::ostringstream error_message;
56  error_message << "Strain matrix is " << strain.ncol() << " x "
57  << strain.nrow() << ", but dimension of the equations is "
58  << DIM << std::endl;
59  throw OomphLibError(
60  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
61  }
62 
63  // Find out how many position types there are
64  unsigned n_position_type = this->nnodal_position_type();
65 
66  if (n_position_type != 1)
67  {
68  throw OomphLibError("TimeHarmonicLinearElasticity is not yet implemented "
69  "for more than one position type",
70  OOMPH_CURRENT_FUNCTION,
71  OOMPH_EXCEPTION_LOCATION);
72  }
73 #endif
74 
75 
76  // Find out how many nodes there are in the element
77  unsigned n_node = nnode();
78 
79  // Find the indices at which the local velocities are stored
80  std::complex<unsigned> u_nodal_index[DIM];
81  for (unsigned i = 0; i < DIM; i++)
82  {
83  u_nodal_index[i] = u_index_time_harmonic_linear_elasticity(i);
84  }
85 
86  // Set up memory for the shape and derivative functions
87  Shape psi(n_node);
88  DShape dpsidx(n_node, DIM);
89 
90  // Call the derivatives of the shape functions
91  (void)dshape_eulerian(s, psi, dpsidx);
92 
93  // Calculate interpolated values of the derivative of global position
94  DenseMatrix<std::complex<double>> interpolated_dudx(
95  DIM, DIM, std::complex<double>(0.0, 0.0));
96 
97  // Loop over nodes
98  for (unsigned l = 0; l < n_node; l++)
99  {
100  // Loop over velocity components
101  for (unsigned i = 0; i < DIM; i++)
102  {
103  // Get the nodal value
104  const std::complex<double> u_value =
105  std::complex<double>(this->nodal_value(l, u_nodal_index[i].real()),
106  this->nodal_value(l, u_nodal_index[i].imag()));
107 
108  // Loop over derivative directions
109  for (unsigned j = 0; j < DIM; j++)
110  {
111  interpolated_dudx(i, j) += u_value * dpsidx(l, j);
112  }
113  }
114  }
115 
116  /// Now fill in the entries of the strain tensor
117  for (unsigned i = 0; i < DIM; i++)
118  {
119  // Do upper half of matrix
120  // Note that j must be signed here for the comparison test to work
121  // Also i must be cast to an int
122  for (int j = (DIM - 1); j >= static_cast<int>(i); j--)
123  {
124  // Off diagonal terms
125  if (static_cast<int>(i) != j)
126  {
127  strain(i, j) =
128  0.5 * (interpolated_dudx(i, j) + interpolated_dudx(j, i));
129  }
130  // Diagonal terms will including growth factor when it comes back in
131  else
132  {
133  strain(i, i) = interpolated_dudx(i, i);
134  }
135  }
136  // Matrix is symmetric so just copy lower half
137  for (int j = (i - 1); j >= 0; j--)
138  {
139  strain(i, j) = strain(j, i);
140  }
141  }
142  }
143 
144 
145  //======================================================================
146  /// Compute the Cauchy stress tensor at local coordinate s for
147  /// displacement formulation.
148  //======================================================================
149  template<unsigned DIM>
151  const Vector<double>& s, DenseMatrix<std::complex<double>>& stress) const
152  {
153 #ifdef PARANOID
154  if ((stress.ncol() != DIM) || (stress.nrow() != DIM))
155  {
156  std::ostringstream error_message;
157  error_message << "Stress matrix is " << stress.ncol() << " x "
158  << stress.nrow() << ", but dimension of the equations is "
159  << DIM << std::endl;
160  throw OomphLibError(
161  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
162  }
163 #endif
164 
165  // Get strain
166  DenseMatrix<std::complex<double>> strain(DIM, DIM);
167  this->get_strain(s, strain);
168 
169  // Now fill in the entries of the stress tensor without exploiting
170  // symmetry -- sorry too lazy. This fct is only used for
171  // postprocessing anyway...
172  for (unsigned i = 0; i < DIM; i++)
173  {
174  for (unsigned j = 0; j < DIM; j++)
175  {
176  stress(i, j) = 0.0;
177  for (unsigned k = 0; k < DIM; k++)
178  {
179  for (unsigned l = 0; l < DIM; l++)
180  {
181  stress(i, j) += this->E(i, j, k, l) * strain(k, l);
182  }
183  }
184  }
185  }
186  }
187 
188 
189  //=======================================================================
190  /// Compute the residuals for the linear elasticity equations in
191  /// cartesian coordinates. Flag indicates if we want Jacobian too.
192  //=======================================================================
193  template<unsigned DIM>
196  Vector<double>& residuals, DenseMatrix<double>& jacobian, unsigned flag)
197  {
198  // Find out how many nodes there are
199  unsigned n_node = this->nnode();
200 
201 #ifdef PARANOID
202  // Find out how many positional dofs there are
203  unsigned n_position_type = this->nnodal_position_type();
204 
205  if (n_position_type != 1)
206  {
207  throw OomphLibError("TimeHarmonicLinearElasticity is not yet implemented "
208  "for more than one position type",
209  OOMPH_CURRENT_FUNCTION,
210  OOMPH_EXCEPTION_LOCATION);
211  }
212 
213  // Throw and error if an elasticity tensor has not been set
214  if (this->Elasticity_tensor_pt == 0)
215  {
216  throw OomphLibError("No elasticity tensor set.",
217  OOMPH_CURRENT_FUNCTION,
218  OOMPH_EXCEPTION_LOCATION);
219  }
220 #endif
221 
222  // Find the indices at which the local velocities are stored
223  std::complex<unsigned> u_nodal_index[DIM];
224  for (unsigned i = 0; i < DIM; i++)
225  {
226  u_nodal_index[i] = this->u_index_time_harmonic_linear_elasticity(i);
227  }
228 
229 
230  // Square of non-dimensional frequency
231  const double omega_sq_local = this->omega_sq();
232 
233  // Set up memory for the shape functions
234  Shape psi(n_node);
235  DShape dpsidx(n_node, DIM);
236 
237  // Set the value of Nintpt -- the number of integration points
238  unsigned n_intpt = this->integral_pt()->nweight();
239 
240  // Set the vector to hold the local coordinates in the element
241  Vector<double> s(DIM);
242 
243  // Integer to store the local equation number
244  int local_eqn = 0, local_unknown = 0;
245 
246  // Loop over the integration points
247  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
248  {
249  // Assign the values of s
250  for (unsigned i = 0; i < DIM; ++i)
251  {
252  s[i] = this->integral_pt()->knot(ipt, i);
253  }
254 
255  // Get the integral weight
256  double w = this->integral_pt()->weight(ipt);
257 
258  // Call the derivatives of the shape functions (and get Jacobian)
259  double J = this->dshape_eulerian_at_knot(ipt, psi, dpsidx);
260 
261  // Storage for Eulerian coordinates (initialised to zero)
262  Vector<double> interpolated_x(DIM, 0.0);
263 
264  // Displacement
265  Vector<std::complex<double>> interpolated_u(
266  DIM, std::complex<double>(0.0, 0.0));
267 
268  // Calculate interpolated values of the derivative of global position
269  // wrt lagrangian coordinates
270  DenseMatrix<std::complex<double>> interpolated_dudx(
271  DIM, DIM, std::complex<double>(0.0, 0.0));
272 
273  // Calculate displacements and derivatives
274  for (unsigned l = 0; l < n_node; l++)
275  {
276  // Loop over displacement components (deformed position)
277  for (unsigned i = 0; i < DIM; i++)
278  {
279  // Calculate the Lagrangian coordinates and the accelerations
280  interpolated_x[i] += this->raw_nodal_position(l, i) * psi(l);
281 
282 
283  // Get the nodal displacements
284  const std::complex<double> u_value = std::complex<double>(
285  this->raw_nodal_value(l, u_nodal_index[i].real()),
286  this->raw_nodal_value(l, u_nodal_index[i].imag()));
287 
288  interpolated_u[i] += u_value * psi(l);
289 
290  // Loop over derivative directions
291  for (unsigned j = 0; j < DIM; j++)
292  {
293  interpolated_dudx(i, j) += u_value * dpsidx(l, j);
294  }
295  }
296  }
297 
298  // Get body force at current time
300  this->body_force(interpolated_x, b);
301 
302  // Premultiply the weights and the Jacobian
303  double W = w * J;
304 
305  //=====EQUATIONS OF LINEAR ELASTICITY ========
306 
307  // Loop over the test functions, nodes of the element
308  for (unsigned l = 0; l < n_node; l++)
309  {
310  // Loop over the displacement components
311  for (unsigned a = 0; a < DIM; a++)
312  {
313  // Get the REAL equation number
314  local_eqn = this->nodal_local_eqn(l, u_nodal_index[a].real());
315 
316  /*IF it's not a boundary condition*/
317  if (local_eqn >= 0)
318  {
319  // Acceleration and body force
320  residuals[local_eqn] +=
321  (-omega_sq_local * interpolated_u[a].real() - b[a].real()) *
322  psi(l) * W;
323 
324  // Stress term
325  for (unsigned b = 0; b < DIM; b++)
326  {
327  for (unsigned c = 0; c < DIM; c++)
328  {
329  for (unsigned d = 0; d < DIM; d++)
330  {
331  // Add the stress terms to the residuals
332  residuals[local_eqn] += this->E(a, b, c, d) *
333  interpolated_dudx(c, d).real() *
334  dpsidx(l, b) * W;
335  }
336  }
337  }
338 
339  // Jacobian entries
340  if (flag)
341  {
342  // Loop over the displacement basis functions again
343  for (unsigned l2 = 0; l2 < n_node; l2++)
344  {
345  // Loop over the displacement components again
346  for (unsigned c = 0; c < DIM; c++)
347  {
348  local_unknown =
349  this->nodal_local_eqn(l2, u_nodal_index[c].real());
350  // If it's not pinned
351  if (local_unknown >= 0)
352  {
353  // Inertial term
354  if (a == c)
355  {
356  jacobian(local_eqn, local_unknown) -=
357  omega_sq_local * psi(l) * psi(l2) * W;
358  }
359 
360  // Stress term
361  for (unsigned b = 0; b < DIM; b++)
362  {
363  for (unsigned d = 0; d < DIM; d++)
364  {
365  // Add the contribution to the Jacobian matrix
366  jacobian(local_eqn, local_unknown) +=
367  this->E(a, b, c, d) * dpsidx(l2, d) * dpsidx(l, b) *
368  W;
369  }
370  }
371  } // End of if not boundary condition
372  }
373  }
374  } // End of jacobian calculation
375 
376  } // End of if not boundary condition for real eqn
377 
378 
379  // Get the IMAG equation number
380  local_eqn = this->nodal_local_eqn(l, u_nodal_index[a].imag());
381 
382  /*IF it's not a boundary condition*/
383  if (local_eqn >= 0)
384  {
385  // Acceleration and body force
386  residuals[local_eqn] +=
387  (-omega_sq_local * interpolated_u[a].imag() - b[a].imag()) *
388  psi(l) * W;
389 
390  // Stress term
391  for (unsigned b = 0; b < DIM; b++)
392  {
393  for (unsigned c = 0; c < DIM; c++)
394  {
395  for (unsigned d = 0; d < DIM; d++)
396  {
397  // Add the stress terms to the residuals
398  residuals[local_eqn] += this->E(a, b, c, d) *
399  interpolated_dudx(c, d).imag() *
400  dpsidx(l, b) * W;
401  }
402  }
403  }
404 
405  // Jacobian entries
406  if (flag)
407  {
408  // Loop over the displacement basis functions again
409  for (unsigned l2 = 0; l2 < n_node; l2++)
410  {
411  // Loop over the displacement components again
412  for (unsigned c = 0; c < DIM; c++)
413  {
414  local_unknown =
415  this->nodal_local_eqn(l2, u_nodal_index[c].imag());
416  // If it's not pinned
417  if (local_unknown >= 0)
418  {
419  // Inertial term
420  if (a == c)
421  {
422  jacobian(local_eqn, local_unknown) -=
423  omega_sq_local * psi(l) * psi(l2) * W;
424  }
425 
426  // Stress term
427  for (unsigned b = 0; b < DIM; b++)
428  {
429  for (unsigned d = 0; d < DIM; d++)
430  {
431  // Add the contribution to the Jacobian matrix
432  jacobian(local_eqn, local_unknown) +=
433  this->E(a, b, c, d) * dpsidx(l2, d) * dpsidx(l, b) *
434  W;
435  }
436  }
437  } // End of if not boundary condition
438  }
439  }
440  } // End of jacobian calculation
441 
442  } // End of if not boundary condition for imag eqn
443 
444  } // End of loop over coordinate directions
445  } // End of loop over shape functions
446  } // End of loop over integration points
447  }
448 
449  //=======================================================================
450  /// Output exact solution x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i]
451  //=======================================================================
452  template<unsigned DIM>
454  std::ostream& outfile,
455  const unsigned& nplot,
457  {
458  // Vector of local coordinates
459  Vector<double> s(DIM);
460 
461  // Vector for coordintes
462  Vector<double> x(DIM);
463 
464  // Tecplot header info
465  outfile << this->tecplot_zone_string(nplot);
466 
467  // Exact solution Vector
468  Vector<double> exact_soln(2 * DIM);
469 
470  // Loop over plot points
471  unsigned num_plot_points = this->nplot_points(nplot);
472  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
473  {
474  // Get local coordinates of plot point
475  this->get_s_plot(iplot, nplot, s);
476 
477  // Get x position as Vector
478  this->interpolated_x(s, x);
479 
480  // Get exact solution at this point
481  (*exact_soln_pt)(x, exact_soln);
482 
483  // Output x,y,...,u_exact,...
484  for (unsigned i = 0; i < DIM; i++)
485  {
486  outfile << x[i] << " ";
487  }
488  for (unsigned i = 0; i < 2 * DIM; i++)
489  {
490  outfile << exact_soln[i] << " ";
491  }
492  outfile << std::endl;
493  }
494 
495  // Write tecplot footer (e.g. FE connectivity lists)
496  this->write_tecplot_zone_footer(outfile, nplot);
497  }
498 
499 
500  //=======================================================================
501  /// Output: x,y,[z],u,v,[w]
502  //=======================================================================
503  template<unsigned DIM>
505  const unsigned& nplot)
506  {
507  // Set output Vector
508  Vector<double> s(DIM);
509  Vector<double> x(DIM);
511 
512  // Tecplot header info
513  outfile << this->tecplot_zone_string(nplot);
514 
515  // Loop over plot points
516  unsigned num_plot_points = this->nplot_points(nplot);
517  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
518  {
519  // Get local coordinates of plot point
520  this->get_s_plot(iplot, nplot, s);
521 
522  // Get Eulerian coordinates and displacements
523  this->interpolated_x(s, x);
524  this->interpolated_u_time_harmonic_linear_elasticity(s, u);
525 
526  // Output the x,y,..
527  for (unsigned i = 0; i < DIM; i++)
528  {
529  outfile << x[i] << " ";
530  }
531 
532  // Output u,v,..
533  for (unsigned i = 0; i < DIM; i++)
534  {
535  outfile << u[i].real() << " ";
536  }
537 
538  // Output u,v,..
539  for (unsigned i = 0; i < DIM; i++)
540  {
541  outfile << u[i].imag() << " ";
542  }
543 
544  outfile << std::endl;
545  }
546 
547  // Write tecplot footer (e.g. FE connectivity lists)
548  this->write_tecplot_zone_footer(outfile, nplot);
549  }
550 
551 
552  //=======================================================================
553  /// C-style output: x,y,[z],u,v,[w]
554  //=======================================================================
555  template<unsigned DIM>
557  const unsigned& nplot)
558  {
559  // Vector of local coordinates
560  Vector<double> s(DIM);
561 
562  // Tecplot header info
563  fprintf(file_pt, "%s", this->tecplot_zone_string(nplot).c_str());
564 
565  // Loop over plot points
566  unsigned num_plot_points = this->nplot_points(nplot);
567  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
568  {
569  // Get local coordinates of plot point
570  this->get_s_plot(iplot, nplot, s);
571 
572  // Coordinates
573  for (unsigned i = 0; i < DIM; i++)
574  {
575  fprintf(file_pt, "%g ", this->interpolated_x(s, i));
576  }
577 
578  // Displacement
579  for (unsigned i = 0; i < DIM; i++)
580  {
581  fprintf(
582  file_pt,
583  "%g ",
584  this->interpolated_u_time_harmonic_linear_elasticity(s, i).real());
585  }
586  for (unsigned i = 0; i < DIM; i++)
587  {
588  fprintf(
589  file_pt,
590  "%g ",
591  this->interpolated_u_time_harmonic_linear_elasticity(s, i).imag());
592  }
593  }
594  fprintf(file_pt, "\n");
595 
596  // Write tecplot footer (e.g. FE connectivity lists)
597  this->write_tecplot_zone_footer(file_pt, nplot);
598  }
599 
600 
601  //=======================================================================
602  /// Compute norm of the solution
603  //=======================================================================
604  template<unsigned DIM>
606  {
607  // Initialise
608  norm = 0.0;
609 
610  // Vector of local coordinates
611  Vector<double> s(DIM);
612 
613  // Vector for coordintes
614  Vector<double> x(DIM);
615 
616  // Displacement vector
617  Vector<std::complex<double>> disp(DIM);
618 
619  // Find out how many nodes there are in the element
620  unsigned n_node = this->nnode();
621 
622  Shape psi(n_node);
623 
624  // Set the value of n_intpt
625  unsigned n_intpt = this->integral_pt()->nweight();
626 
627  // Loop over the integration points
628  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
629  {
630  // Assign values of s
631  for (unsigned i = 0; i < DIM; i++)
632  {
633  s[i] = this->integral_pt()->knot(ipt, i);
634  }
635 
636  // Get the integral weight
637  double w = this->integral_pt()->weight(ipt);
638 
639  // Get jacobian of mapping
640  double J = this->J_eulerian(s);
641 
642  // Premultiply the weights and the Jacobian
643  double W = w * J;
644 
645  // Get FE function value
646  this->interpolated_u_time_harmonic_linear_elasticity(s, disp);
647 
648  // Add to norm
649  for (unsigned ii = 0; ii < DIM; ii++)
650  {
651  norm += (disp[ii].real() * disp[ii].real() +
652  disp[ii].imag() * disp[ii].imag()) *
653  W;
654  }
655  }
656  }
657 
658 
659  // Instantiate the required elements
662 
666 
667 
668 } // namespace oomph
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition: shape.h:278
//////////////////////////////////////////////////////////////////////////// ////////////////////////...
Definition: matrices.h:386
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
An OomphLibError object which should be thrown when an run-time error is encountered....
///////////////////////////////////////////////////////////////////// ///////////////////////////////...
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
A base class for elements that solve the equations of time-harmonic linear elasticity in Cartesian co...
static double Default_omega_sq_value
Static default value for square of frequency.
void get_strain(const Vector< double > &s, DenseMatrix< std::complex< double >> &strain) const
Return the strain tensor.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...
virtual void fill_in_generic_contribution_to_residuals_time_harmonic_linear_elasticity(Vector< double > &residuals, DenseMatrix< double > &jacobian, unsigned flag)
Private helper function to compute residuals and (if requested via flag) also the Jacobian matrix.
void compute_norm(double &norm)
Compute norm of solution: square of the L2 norm.
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
void output(std::ostream &outfile)
Output: x,y,[z],u_r,v_r,[w_r],u_i,v_i,[w_i].
void get_stress(const Vector< double > &s, DenseMatrix< std::complex< double >> &sigma) const
Return the Cauchy stress tensor, as calculated from the elasticity tensor at specified local coordina...
//////////////////////////////////////////////////////////////////// ////////////////////////////////...