axisym_fvk_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 axisym FoepplvonKarman elements
27 
28 #include "axisym_fvk_elements.h"
29 
30 namespace oomph
31 {
32  //======================================================================
33  /// Set the data for the number of Variables at each node: 6
34  //======================================================================
35  template<unsigned NNODE_1D>
37 
38 
39  //======================================================================
40  /// Default value physical constants
41  //======================================================================
43 
44 
45  //======================================================================
46  /// Compute contribution to element residual Vector
47  ///
48  /// Pure version without hanging nodes
49  //======================================================================
51  Vector<double>& residuals)
52  {
53  // Find out how many nodes there are
54  const unsigned n_node = nnode();
55 
56  // Set up memory for the shape and test functions
57  Shape psi(n_node), test(n_node);
58  DShape dpsidr(n_node, 1), dtestdr(n_node, 1);
59 
60  // Indices at which the unknowns are stored
61  const unsigned w_nodal_index = nodal_index_fvk(0);
62  const unsigned laplacian_w_nodal_index = nodal_index_fvk(1);
63  const unsigned phi_nodal_index = nodal_index_fvk(2);
64  const unsigned laplacian_phi_nodal_index = nodal_index_fvk(3);
65  const unsigned smooth_dwdr_nodal_index = nodal_index_fvk(4);
66  const unsigned smooth_dphidr_nodal_index = nodal_index_fvk(5);
67 
68  // Set the value of n_intpt
69  const unsigned n_intpt = integral_pt()->nweight();
70 
71  // Integers to store the local equation numbers
72  int local_eqn = 0;
73 
74  // Loop over the integration points
75  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
76  {
77  // Get the integral weight
78  double w = integral_pt()->weight(ipt);
79 
80  // Call the derivatives of the shape and test functions
82  ipt, psi, dpsidr, test, dtestdr);
83 
84  // Allocate and initialise to zero storage for the interpolated values
85  double interpolated_r = 0;
86 
87  double interpolated_w = 0;
88  double interpolated_laplacian_w = 0;
89  double interpolated_phi = 0;
90  double interpolated_laplacian_phi = 0;
91 
92  double interpolated_dwdr = 0;
93  double interpolated_dlaplacian_wdr = 0;
94  double interpolated_dphidr = 0;
95  double interpolated_dlaplacian_phidr = 0;
96 
97  double interpolated_smooth_dwdr = 0;
98  double interpolated_smooth_dphidr = 0;
99  double interpolated_continuous_d2wdr2 = 0;
100  double interpolated_continuous_d2phidr2 = 0;
101 
102  // Calculate function values and derivatives:
103  //-----------------------------------------
104  Vector<double> nodal_value(6, 0.0);
105  // Loop over nodes
106  for (unsigned l = 0; l < n_node; l++)
107  {
108  // Get the nodal values
109  nodal_value[0] = raw_nodal_value(l, w_nodal_index);
110  nodal_value[1] = raw_nodal_value(l, laplacian_w_nodal_index);
111 
113  {
114  nodal_value[2] = raw_nodal_value(l, phi_nodal_index);
115  nodal_value[3] = raw_nodal_value(l, laplacian_phi_nodal_index);
116  nodal_value[4] = raw_nodal_value(l, smooth_dwdr_nodal_index);
117  nodal_value[5] = raw_nodal_value(l, smooth_dphidr_nodal_index);
118  }
119 
120  // Add contributions from current node/shape function
121  interpolated_w += nodal_value[0] * psi(l);
122  interpolated_laplacian_w += nodal_value[1] * psi(l);
123 
125  {
126  interpolated_phi += nodal_value[2] * psi(l);
127  interpolated_laplacian_phi += nodal_value[3] * psi(l);
128  interpolated_smooth_dwdr += nodal_value[4] * psi(l);
129  interpolated_smooth_dphidr += nodal_value[5] * psi(l);
130 
131  interpolated_continuous_d2wdr2 += nodal_value[4] * dpsidr(l, 0);
132  interpolated_continuous_d2phidr2 += nodal_value[5] * dpsidr(l, 0);
133  }
134 
135  interpolated_r += raw_nodal_position(l, 0) * psi(l);
136  interpolated_dwdr += nodal_value[0] * dpsidr(l, 0);
137  interpolated_dlaplacian_wdr += nodal_value[1] * dpsidr(l, 0);
138 
140  {
141  interpolated_dphidr += nodal_value[2] * dpsidr(l, 0);
142  interpolated_dlaplacian_phidr += nodal_value[3] * dpsidr(l, 0);
143  }
144  } // End of loop over the nodes
145 
146 
147  // Premultiply the weights and the Jacobian
148  double W = w * interpolated_r * J;
149 
150  // Get pressure function
151  //-------------------
152  double pressure;
153  get_pressure_fvk(ipt, interpolated_r, pressure);
154 
155  double airy_forcing;
156  get_airy_forcing_fvk(ipt, interpolated_r, airy_forcing);
157 
158 
159  // Assemble residuals and Jacobian
160  //--------------------------------
161 
162  // Loop over the test functions
163  for (unsigned l = 0; l < n_node; l++)
164  {
165  // Get the local equation
166  local_eqn = nodal_local_eqn(l, w_nodal_index);
167 
168  // IF it's not a boundary condition
169  if (local_eqn >= 0)
170  {
171  residuals[local_eqn] += pressure * test(l) * W;
172 
173  // Reduced order biharmonic operator
174  residuals[local_eqn] +=
175  interpolated_dlaplacian_wdr * dtestdr(l, 0) * W;
176 
178  {
179  // Monge-Ampere part
180  residuals[local_eqn] +=
181  eta() *
182  (interpolated_continuous_d2wdr2 * interpolated_dphidr +
183  interpolated_dwdr * interpolated_continuous_d2phidr2) /
184  interpolated_r * test(l) * W;
185  }
186  }
187 
188  // Get the local equation
189  local_eqn = nodal_local_eqn(l, laplacian_w_nodal_index);
190 
191  // IF it's not a boundary condition
192  if (local_eqn >= 0)
193  {
194  // The coupled Poisson equations for the biharmonic operator
195  residuals[local_eqn] += interpolated_laplacian_w * test(l) * W;
196  residuals[local_eqn] += interpolated_dwdr * dtestdr(l, 0) * W;
197  }
198 
199  // Get the local equation
200  local_eqn = nodal_local_eqn(l, phi_nodal_index);
201 
202  // IF it's not a boundary condition
203  if (local_eqn >= 0)
204  {
205  residuals[local_eqn] += airy_forcing * test(l) * W;
206 
207  // Reduced order biharmonic operator
208  residuals[local_eqn] +=
209  interpolated_dlaplacian_phidr * dtestdr(l, 0) * W;
210 
211  // Monge-Ampere part
212  residuals[local_eqn] -= interpolated_continuous_d2wdr2 *
213  interpolated_dwdr / interpolated_r * test(l) *
214  W;
215  }
216 
217  // Get the local equation
218  local_eqn = nodal_local_eqn(l, laplacian_phi_nodal_index);
219 
220  // IF it's not a boundary condition
221  if (local_eqn >= 0)
222  {
223  // The coupled Poisson equations for the biharmonic operator
224  residuals[local_eqn] += interpolated_laplacian_phi * test(l) * W;
225  residuals[local_eqn] += interpolated_dphidr * dtestdr(l, 0) * W;
226  }
227 
228  // Residuals for the smooth derivatives
229  local_eqn = nodal_local_eqn(l, smooth_dwdr_nodal_index);
230  if (local_eqn >= 0)
231  {
232  residuals[local_eqn] +=
233  (interpolated_dwdr - interpolated_smooth_dwdr) * test(l) * W;
234  }
235 
236  local_eqn = nodal_local_eqn(l, smooth_dphidr_nodal_index);
237  if (local_eqn >= 0)
238  {
239  residuals[local_eqn] +=
240  (interpolated_dphidr - interpolated_smooth_dphidr) * test(l) * W;
241  }
242 
243  } // End of loop over test functions
244  } // End of loop over integration points
245  }
246 
247 
248  //======================================================================
249  /// Self-test: Return 0 for OK
250  //======================================================================
252  {
253  bool passed = true;
254 
255  // Check lower-level stuff
256  if (FiniteElement::self_test() != 0)
257  {
258  passed = false;
259  }
260 
261  // Return verdict
262  if (passed)
263  {
264  return 0;
265  }
266  else
267  {
268  return 1;
269  }
270  }
271 
272 
273  //======================================================================
274  /// Compute in-plane stresses.
275  //======================================================================
277  const Vector<double>& s, double& sigma_r_r, double& sigma_phi_phi) const
278  {
279  // No in plane stresses if linear bending
281  {
282  sigma_r_r = 0.0;
283  sigma_phi_phi = 0.0;
284  }
285  else
286  {
287  // Get shape fcts and derivs
288  unsigned n_dim = this->dim();
289  unsigned n_node = this->nnode();
290  Shape psi(n_node);
291  DShape dpsi_dr(n_node, n_dim);
292 
293  // Get shape fcts and derivs
294  dshape_eulerian(s, psi, dpsi_dr);
295 
296  // Indices at which the unknowns are stored
297  const unsigned smooth_dphi_dr_nodal_index = nodal_index_fvk(5);
298 
299  // Allocate and initialise to zero storage for the interpolated values
300  double interpolated_r = 0;
301  double interpolated_dphi_dr = 0;
302  double interpolated_continuous_d2phi_dr2 = 0;
303 
304 
305  // Calculate function values and derivatives:
306  //-----------------------------------------
307  // Loop over nodes
308  for (unsigned l = 0; l < n_node; l++)
309  {
310  // Add contributions from current node/shape function
311  interpolated_r += raw_nodal_position(l, 0) * psi(l);
312  interpolated_dphi_dr +=
313  this->nodal_value(l, smooth_dphi_dr_nodal_index) * psi(l);
314  interpolated_continuous_d2phi_dr2 +=
315  this->nodal_value(l, smooth_dphi_dr_nodal_index) * dpsi_dr(l, 0);
316  } // End of loop over nodes
317 
318  // Compute stresses
319  if (interpolated_r == 0.0)
320  {
321  // To avoid dividing by zero we take advantage of axisymmetry at the
322  // origin: l'Hopitals rule => dphi_dr/r -> d2phi_dr2
323  sigma_r_r = interpolated_continuous_d2phi_dr2;
324  }
325  else
326  {
327  sigma_r_r = interpolated_dphi_dr / interpolated_r;
328  }
329  sigma_phi_phi = interpolated_continuous_d2phi_dr2;
330 
331  } // End if
332 
333  } // End of interpolated_stress function
334 
335 
336  //======================================================================
337  /// Output function:
338  /// r, w, sigma_r_r, sigma_phi_phi
339  /// nplot points
340  //======================================================================
341  void AxisymFoepplvonKarmanEquations::output(std::ostream& outfile,
342  const unsigned& nplot)
343  {
344  // Vector of local coordinates
345  Vector<double> s(1);
346 
347  // Tecplot header info
348  outfile << "ZONE\n";
349 
350  // Loop over plot points
351  unsigned num_plot_points = nplot_points(nplot);
352  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
353  {
354  // Get local coordinates of plot point
355  get_s_plot(iplot, nplot, s);
356 
357  // Get stress
358  double sigma_r_r = 0.0;
359  double sigma_phi_phi = 0.0;
360  interpolated_stress(s, sigma_r_r, sigma_phi_phi);
361 
362  // Output interpolated global position, deflection and stress
363  outfile << interpolated_x(s, 0) << " " << interpolated_w_fvk(s) << " "
364  << sigma_r_r << " " << sigma_phi_phi << std::endl;
365  }
366  }
367 
368 
369  //======================================================================
370  /// C-style output function:
371  /// r,w
372  /// nplot points
373  //======================================================================
374  void AxisymFoepplvonKarmanEquations::output(FILE* file_pt,
375  const unsigned& nplot)
376  {
377  // Vector of local coordinates
378  Vector<double> s(1);
379 
380  // Tecplot header info
381  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
382 
383  // Loop over plot points
384  unsigned num_plot_points = nplot_points(nplot);
385  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
386  {
387  // Get local coordinates of plot point
388  get_s_plot(iplot, nplot, s);
389 
390  fprintf(file_pt, "%g ", interpolated_x(s, 0));
391  fprintf(file_pt, "%g \n", interpolated_w_fvk(s));
392  }
393 
394  // Write tecplot footer (e.g. FE connectivity lists)
395  write_tecplot_zone_footer(file_pt, nplot);
396  }
397 
398 
399  //======================================================================
400  /// Output exact solution
401  ///
402  /// Solution is provided via function pointer.
403  /// Plot at a given number of plot points.
404  ///
405  /// r,w_exact
406  //======================================================================
408  std::ostream& outfile,
409  const unsigned& nplot,
411  {
412  // Vector of local coordinates
413  Vector<double> s(1);
414 
415  // Vector for coordinates
416  Vector<double> r(1);
417 
418  // Tecplot header info
419  outfile << tecplot_zone_string(nplot);
420 
421  // Exact solution Vector (here a scalar)
422  // Vector<double> exact_soln(1);
423  Vector<double> exact_soln(1);
424 
425  // Loop over plot points
426  unsigned num_plot_points = nplot_points(nplot);
427  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
428  {
429  // Get local coordinates of plot point
430  get_s_plot(iplot, nplot, s);
431 
432  // Get r position as Vector
433  interpolated_x(s, r);
434 
435  // Get exact solution at this point
436  (*exact_soln_pt)(r, exact_soln);
437 
438  // Output r,w_exact
439  outfile << r[0] << " ";
440  outfile << exact_soln[0] << std::endl;
441  }
442 
443  // Write tecplot footer (e.g. FE connectivity lists)
444  write_tecplot_zone_footer(outfile, nplot);
445  }
446 
447 
448  //======================================================================
449  /// Validate against exact solution
450  ///
451  /// Solution is provided via function pointer.
452  /// Plot error at a given number of plot points.
453  ///
454  //======================================================================
456  std::ostream& outfile,
458  double& error,
459  double& norm)
460  {
461  // Initialise
462  error = 0.0;
463  norm = 0.0;
464 
465  // Vector of local coordinates
466  Vector<double> s(1);
467 
468  // Vector for coordintes
469  Vector<double> r(1);
470 
471  // Find out how many nodes there are in the element
472  unsigned n_node = nnode();
473 
474  Shape psi(n_node);
475 
476  // Set the value of n_intpt
477  unsigned n_intpt = integral_pt()->nweight();
478 
479  // Tecplot
480  outfile << "ZONE" << std::endl;
481 
482  // Exact solution Vector (here a scalar)
483  // Vector<double> exact_soln(1);
484  Vector<double> exact_soln(1);
485 
486  // Loop over the integration points
487  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
488  {
489  // Assign values of s
490  s[0] = integral_pt()->knot(ipt, 0);
491 
492  // Get the integral weight
493  double w = integral_pt()->weight(ipt);
494 
495  // Get jacobian of mapping
496  double J = J_eulerian(s);
497 
498  // Premultiply the weights and the Jacobian
499  double W = w * J;
500 
501  // Get r position as Vector
502  interpolated_x(s, r);
503 
504  // Get FE function value
505  double w_fe = interpolated_w_fvk(s);
506 
507  // Get exact solution at this point
508  (*exact_soln_pt)(r, exact_soln);
509 
510  // Output r,error
511  outfile << r[0] << " ";
512  outfile << exact_soln[0] << " " << exact_soln[0] - w_fe << std::endl;
513 
514  // Add to error and norm
515  norm += exact_soln[0] * exact_soln[0] * W;
516  error += (exact_soln[0] - w_fe) * (exact_soln[0] - w_fe) * W;
517  }
518  }
519 
520 
521  //====================================================================
522  // Force build of templates
523  //====================================================================
524  template class AxisymFoepplvonKarmanElement<2>;
525  template class AxisymFoepplvonKarmanElement<3>;
526  template class AxisymFoepplvonKarmanElement<4>;
527 
528 } // namespace oomph
static char t char * s
Definition: cfortran.h:568
static const unsigned Initial_Nvalue
Static int that holds the number of variables at nodes: always the same.
static double Default_Physical_Constant_Value
Default value for physical constants.
void output(std::ostream &outfile)
Output with default number of plot points.
void fill_in_contribution_to_residuals(Vector< double > &residuals)
Fill in the residuals with this element's contribution.
unsigned self_test()
Self-test: Return 0 for OK.
virtual double dshape_and_dtest_eulerian_at_knot_axisym_fvk(const unsigned &ipt, Shape &psi, DShape &dpsidr, Shape &test, DShape &dtestdr) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
virtual unsigned nodal_index_fvk(const unsigned &i=0) const
Return the index at which the i-th unknown value is stored. The default value, i, is appropriate for ...
void interpolated_stress(const Vector< double > &s, double &sigma_r_r, double &sigma_phi_phi) const
Compute in-plane stresses.
bool Linear_bending_model
Flag which stores whether we are using a linear, pure bending model instead of the full non-linear Fo...
const double & eta() const
FvK parameter.
void output_fct(std::ostream &outfile, const unsigned &n_plot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,w_exact at n_plot plot points.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
double interpolated_w_fvk(const Vector< double > &s) const
Return FE representation of transverse displacement.
virtual void get_pressure_fvk(const unsigned &ipt, const double &r, double &pressure) const
Get pressure term at (Eulerian) position r. This function is virtual to allow overloading in multi-ph...
virtual void get_airy_forcing_fvk(const unsigned &ipt, const double &r, double &airy_forcing) const
Get Airy forcing term at (Eulerian) position r. This function is virtual to allow overloading in mult...
virtual double J_eulerian(const Vector< double > &s) const
Return the Jacobian of mapping from local to global coordinates at local position s.
Definition: elements.cc:4133
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition: elements.h:2597
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
Definition: elements.h:3165
virtual double interpolated_x(const Vector< double > &s, const unsigned &i) const
Return FE interpolated coordinate x[i] at local coordinate s.
Definition: elements.cc:3992
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition: elements.h:1436
unsigned dim() const
Return the spatial dimension of the element, i.e. the number of local coordinates required to paramet...
Definition: elements.h:2615
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1763
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1967
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
Definition: elements.h:3152
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
Definition: elements.h:3190
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition: elements.cc:3328
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
Definition: elements.h:2580
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. Do not use the hanging node representation....
Definition: elements.cc:1714
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
Definition: elements.h:3178
virtual unsigned self_test()
Self-test: Check inversion of element & do self-test for GeneralisedElement. Return 0 if OK.
Definition: elements.cc:4470
virtual double knot(const unsigned &i, const unsigned &j) const =0
Return local coordinate s[j] of i-th integration point.
virtual unsigned nweight() const =0
Return the number of integration points of the scheme.
virtual double weight(const unsigned &i) const =0
Return weight of i-th integration point.
//////////////////////////////////////////////////////////////////// ////////////////////////////////...