axisym_displ_based_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-2025 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
29
30
31namespace oomph
32{
33 //======================================================================
34 /// Set the data for the number of Variables at each node - 3
35 //======================================================================
36 template<unsigned NNODE_1D>
38
39
40 //======================================================================
41 /// Compute contribution to element residual Vector
42 ///
43 /// Pure version without hanging nodes
44 //======================================================================
47 {
48 // Find out how many nodes there are
49 const unsigned n_node = nnode();
50
51 // Set up memory for the shape and test functions
54
55 // Set the value of n_intpt
56 const unsigned n_intpt = integral_pt()->nweight();
57
58 // Indices at which the unknowns are stored
59 const unsigned w_nodal_index = nodal_index_fvk(0);
60 const unsigned laplacian_w_nodal_index = nodal_index_fvk(1);
61 const unsigned u_r_nodal_index = nodal_index_fvk(2);
62
63 // Local copy of parameters
64 const double nu_local = nu();
65 const double eta_local = eta();
66
67 // Integers to store the local equation numbers
68 int local_eqn = 0;
69
70 // Loop over the integration points
71 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
72 {
73 // Get the integral weight
74 double w = integral_pt()->weight(ipt);
75
76 // Call the derivatives of the shape and test functions
79
80 // Allocate and initialise to zero storage for the interpolated values
81 double interpolated_r = 0.0;
82
83 double interpolated_w = 0.0;
84 double interpolated_laplacian_w = 0.0;
85 double interpolated_u_r = 0.0;
86
87 double interpolated_dwdr = 0.0;
88 double interpolated_dlaplacian_wdr = 0.0;
89 double interpolated_du_rdr = 0.0;
90
92
93 // Calculate function values and derivatives:
94 //-----------------------------------------
95 // Loop over nodes
96 for (unsigned l = 0; l < n_node; l++)
97 {
98 // Get the nodal values
102
103 // Add contributions from current node/shape function
107
109
113
114 } // End of loop over the nodes
115
116 // Premultiply the weights and the Jacobian
117 double W = w * interpolated_r * J;
118
119 // Get pressure function
120 //---------------------
121 double pressure = 0.0;
123
124 // Determine the stresses
125 //-----------------------
126
127 double sigma_r_r = 0.0;
128 double sigma_phi_phi = 0.0;
129
131 {
132 sigma_r_r =
133 1.0 / (1.0 - nu_local * nu_local) *
136
138 1.0 / (1.0 - nu_local * nu_local) *
142 }
143 else
144 {
145 sigma_r_r = 1.0 / (1.0 - nu_local * nu_local) *
148
149 sigma_phi_phi = 1.0 / (1.0 - nu_local * nu_local) *
152 }
153
154
155 // Assemble residuals and Jacobian:
156 //--------------------------------
157 // Loop over the test functions
158 for (unsigned l = 0; l < n_node; l++)
159 {
160 // Get the local equation
162
163 // IF it's not a boundary condition
164 if (local_eqn >= 0)
165 {
167 (pressure * test(l) +
169 W;
171 {
174 }
175 }
176
177 // Get the local equation
179
180 // IF it's not a boundary condition
181 if (local_eqn >= 0)
182 {
184 (dtestdr(l, 0)) * interpolated_dwdr) *
185 W;
186 }
187
188 // Get the local equation
190
191 // IF it's not a boundary condition
192 if (local_eqn >= 0)
193 {
195 (sigma_r_r * dtestdr(l, 0) +
196 1.0 / interpolated_r * sigma_phi_phi * test(l)) *
197 W;
198 }
199
200 } // End of loop over test functions
201 } // End of loop over integration points
202 }
203
204
205 //======================================================================
206 /// Self-test: Return 0 for OK
207 //======================================================================
209 {
210 bool passed = true;
211
212 // Check lower-level stuff
213 if (FiniteElement::self_test() != 0)
214 {
215 passed = false;
216 }
217
218 // Return verdict
219 if (passed)
220 {
221 return 0;
222 }
223 else
224 {
225 return 1;
226 }
227 }
228
229
230 //======================================================================
231 /// Compute in-plane stresses. Return boolean to indicate success
232 /// (false if attempt to evaluate stresses at zero radius)
233 //======================================================================
235 const Vector<double>& s, double& sigma_r_r, double& sigma_phi_phi) const
236 {
237 // No in plane stresses if linear bending
239 {
240 sigma_r_r = 0.0;
241 sigma_phi_phi = 0.0;
242 }
243 else
244 {
245 // Get shape fcts and derivs
246 unsigned n_dim = this->dim();
247 unsigned n_node = this->nnode();
250
251 // Get shape fcts and derivs
253
254 // Allocate and initialise to zero storage for the interpolated values
255 double interpolated_r = 0.0;
256 double interpolated_u_r = 0.0;
257
258 double interpolated_dwdr = 0.0;
259 double interpolated_du_rdr = 0.0;
260
261 double nu_local = nu();
262
263
264 // Calculate function values and derivatives:
265 //-----------------------------------------
266 // Loop over nodes
267 for (unsigned l = 0; l < n_node; l++)
268 {
269 // Add contributions from current node/shape function
272 this->raw_nodal_value(l, nodal_index_fvk(2)) * psi(l);
274 this->raw_nodal_value(l, nodal_index_fvk(0)) * dpsidr(l, 0);
276 this->raw_nodal_value(l, nodal_index_fvk(2)) * dpsidr(l, 0);
277 } // End of loop over nodes
278
279 // The centre stress must use a different analytical form to avoid
280 // dividing by zero. It can be found from the form below by assuming:
281 // a) u_r=0 at the origin (required under axisymmetry unless there is a
282 // puncture) which leads to the term u_r/r -> du_r/dr as r -> 0 (by
283 // L'Hopital's rule), and
284 // b) given dw/dr=0 at the origin (required for axisymmetry as there can
285 // be no kink in a plate).
286 // Note that this process results in s_{rr}=s_{\phi\phi} at the origin
287 // which is expected as the two coordinate directions are
288 // indistinguishable at the origin under axisymmetry.
289 if (interpolated_r == 0.0)
290 {
291 // Compute the stresses:
294 }
295 else
296 {
297 // Compute the stresses:
298 //---------------------
299 sigma_r_r =
300 1.0 / (1.0 - nu_local * nu_local) *
303
305 1.0 / (1.0 - nu_local * nu_local) *
309 } // End if origin
310 } // End if nonlinear bending
311
312 } // End of interpolated_stress function
313
314 //======================================================================
315 /// Output function:
316 /// r, w, u, sigma_r_r, sigma_phi_phi
317 /// nplot points
318 //======================================================================
320 const unsigned& nplot)
321 {
322 // Vector of local coordinates
323 Vector<double> s(1);
324
325 // Tecplot header info
326 outfile << "ZONE\n";
327
328 // Loop over plot points
330 for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
331 {
332 // Get local coordinates of plot point
334
335 // Get stress
336 double sigma_r_r = 0.0;
337 double sigma_phi_phi = 0.0;
339
340 // Output interpolated global position, displacement and stress
341 outfile << interpolated_x(s, 0) << " " << interpolated_w_fvk(s) << " "
342 << interpolated_u_fvk(s) << " " << sigma_r_r << " "
343 << sigma_phi_phi << std::endl;
344 }
345 }
346
347 //======================================================================
348 /// C-style output function:
349 /// r,w,u
350 /// nplot points
351 //======================================================================
353 const unsigned& nplot)
354 {
355 // Vector of local coordinates
356 Vector<double> s(1);
357
358 // Tecplot header info
360
361 // Loop over plot points
363 for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
364 {
365 // Get local coordinates of plot point
367
368 fprintf(file_pt, "%g ", interpolated_x(s, 0));
371 }
372
373 // Write tecplot footer (e.g. FE connectivity lists)
375 }
376
377
378 //======================================================================
379 /// Output exact solution
380 ///
381 /// Solution is provided via function pointer.
382 /// Plot at a given number of plot points.
383 ///
384 /// r,w_exact
385 //======================================================================
387 std::ostream& outfile,
388 const unsigned& nplot,
390 {
391 // Vector of local coordinates
392 Vector<double> s(1);
393
394 // Vector for coordinates
395 Vector<double> r(1);
396
397 // Tecplot header info
399
400 // Exact solution Vector (here a scalar)
401 // Vector<double> exact_soln(1);
403
404 // Loop over plot points
406 for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
407 {
408 // Get local coordinates of plot point
410
411 // Get r position as Vector
413
414 // Get exact solution at this point
415 (*exact_soln_pt)(r, exact_soln);
416
417 // Output r,w_exact
418 outfile << r[0] << " ";
419 outfile << exact_soln[0] << std::endl;
420 }
421
422 // Write tecplot footer (e.g. FE connectivity lists)
424 }
425
426
427 //======================================================================
428 /// Validate against exact solution
429 ///
430 /// Solution is provided via function pointer.
431 /// Plot error at a given number of plot points.
432 ///
433 //======================================================================
435 std::ostream& outfile,
437 double& error,
438 double& norm)
439 {
440 // Initialise
441 error = 0.0;
442 norm = 0.0;
443
444 // Vector of local coordinates
445 Vector<double> s(1);
446
447 // Vector for coordintes
448 Vector<double> r(1);
449
450 // Find out how many nodes there are in the element
451 unsigned n_node = nnode();
452
454
455 // Set the value of n_intpt
456 unsigned n_intpt = integral_pt()->nweight();
457
458 // Tecplot
459 outfile << "ZONE" << std::endl;
460
461 // Exact solution Vector (here a scalar)
462 // Vector<double> exact_soln(1);
464
465 // Loop over the integration points
466 for (unsigned ipt = 0; ipt < n_intpt; ipt++)
467 {
468 // Assign values of s
469 s[0] = integral_pt()->knot(ipt, 0);
470
471 // Get the integral weight
472 double w = integral_pt()->weight(ipt);
473
474 // Get jacobian of mapping
475 double J = J_eulerian(s);
476
477 // Premultiply the weights and the Jacobian
478 double W = w * J;
479
480 // Get r position as Vector
482
483 // Get FE function value
484 double w_fe = interpolated_w_fvk(s);
485
486 // Get exact solution at this point
487 (*exact_soln_pt)(r, exact_soln);
488
489 // Output r,error
490 outfile << r[0] << " ";
491 outfile << exact_soln[0] << " " << exact_soln[0] - w_fe << std::endl;
492
493 // Add to error and norm
494 norm += exact_soln[0] * exact_soln[0] * W;
495 error += (exact_soln[0] - w_fe) * (exact_soln[0] - w_fe) * W;
496 }
497 }
498
499
500 //====================================================================
501 // Force build of templates
502 //====================================================================
503 template class AxisymFoepplvonKarmanElement<2>;
504 template class AxisymFoepplvonKarmanElement<3>;
505 template class AxisymFoepplvonKarmanElement<4>;
506
507} // 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.
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...
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_u_fvk(const Vector< double > &s) const
Return FE representation of radial displacement.
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...
const double & nu() const
Poisson's ratio.
A Class for the derivatives of shape functions The class design is essentially the same as Shape,...
Definition shape.h:278
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
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition elements.h:1967
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
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.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition shape.h:76
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).