spherical_advection_diffusion_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 Advection Diffusion elements in a spherical
27 // coordinate system
29 
30 namespace oomph
31 {
32  /// 2D Steady Axisymmetric Advection Diffusion elements
33 
34  /// Default value for Peclet number
36 
37 
38  //======================================================================
39  /// Compute element residual Vector and/or element Jacobian matrix
40  ///
41  /// flag=1: compute both
42  /// flag=0: compute only residual Vector
43  ///
44  /// Pure version without hanging nodes
45  //======================================================================
48  Vector<double>& residuals,
49  DenseMatrix<double>& jacobian,
50  DenseMatrix<double>& mass_matrix,
51  unsigned flag)
52  {
53  // Find out how many nodes there are
54  const unsigned n_node = nnode();
55 
56  // Get the nodal index at which the unknown is stored
57  const unsigned u_nodal_index = u_index_spherical_adv_diff();
58 
59  // Set up memory for the shape and test functions
60  Shape psi(n_node), test(n_node);
61  DShape dpsidx(n_node, 2), dtestdx(n_node, 2);
62 
63  // Set the value of n_intpt
64  const unsigned n_intpt = integral_pt()->nweight();
65 
66  // Set the Vector to hold local coordinates
67  Vector<double> s(2);
68 
69  // Get Physical Variables from Element
70  const double scaled_peclet = pe();
71 
72  // Get peclet number multiplied by Strouhal number
73  const double scaled_peclet_st = pe_st();
74 
75  // Integers used to store the local equation number and local unknown
76  // indices for the residuals and jacobians
77  int local_eqn = 0, local_unknown = 0;
78 
79  // Loop over the integration points
80  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
81  {
82  // Assign values of s
83  for (unsigned i = 0; i < 2; i++) s[i] = integral_pt()->knot(ipt, i);
84 
85  // Get the integral weight
86  double w = integral_pt()->weight(ipt);
87 
88  // Call the derivatives of the shape and test functions
90  ipt, psi, dpsidx, test, dtestdx);
91 
92  // Premultiply the weights and the Jacobian
93  double W = w * J;
94 
95  // Calculate local values of the solution and its derivatives
96  // Allocate
97  double interpolated_u = 0.0;
98  double dudt = 0.0;
100  Vector<double> interpolated_dudx(2, 0.0);
101 
102  // Calculate function value and derivatives:
103  //-----------------------------------------
104  // Loop over nodes
105  for (unsigned l = 0; l < n_node; l++)
106  {
107  // Get the value at the node
108  double u_value = raw_nodal_value(l, u_nodal_index);
109  interpolated_u += u_value * psi(l);
110  dudt += du_dt_spherical_adv_diff(l) * psi(l);
111  // Loop over the two coordinates directions
112  for (unsigned j = 0; j < 2; j++)
113  {
114  interpolated_x[j] += raw_nodal_position(l, j) * psi(l);
115  interpolated_dudx[j] += u_value * dpsidx(l, j);
116  }
117  }
118 
119  // Get source function
120  //-------------------
121  double source;
123 
124 
125  // Get wind
126  //--------
127  Vector<double> wind(3);
129 
130  // r is the first position component
131  double r = interpolated_x[0];
132  // theta is the second position component
133  double sin_th = sin(interpolated_x[1]);
134  // dS is the area weighting
135  double dS = r * r * sin_th;
136 
137  // TEMPERATURE EQUATION (Neglected viscous dissipation)
138  // Assemble residuals and Jacobian
139  //--------------------------------
140 
141  // Loop over the test functions
142  for (unsigned l = 0; l < n_node; l++)
143  {
144  // Set the local equation number
145  local_eqn = nodal_local_eqn(l, u_nodal_index);
146 
147  /*IF it's not a boundary condition*/
148  if (local_eqn >= 0)
149  {
150  // Add body force/source term
151  residuals[local_eqn] -=
152  (scaled_peclet_st * dudt + source) * dS * test(l) * W;
153 
154  // The Advection Diffusion bit itself
155  residuals[local_eqn] -=
156  // radial terms
157  (dS * interpolated_dudx[0] *
158  (scaled_peclet * wind[0] * test(l) + dtestdx(l, 0)) +
159  // azimuthal terms
160  (sin_th * interpolated_dudx[1] *
161  (r * scaled_peclet * wind[1] * test(l) + dtestdx(l, 1)))) *
162  W;
163 
164 
165  // Calculate the jacobian
166  //-----------------------
167  if (flag)
168  {
169  // Loop over the velocity shape functions again
170  for (unsigned l2 = 0; l2 < n_node; l2++)
171  {
172  // Set the number of the unknown
173  local_unknown = nodal_local_eqn(l2, u_nodal_index);
174 
175  // If at a non-zero degree of freedom add in the entry
176  if (local_unknown >= 0)
177  {
178  // Mass matrix term
179  jacobian(local_eqn, local_unknown) -=
180  scaled_peclet_st * test(l) * psi(l2) *
181  node_pt(l2)->time_stepper_pt()->weight(1, 0) * dS * W;
182 
183  // add the mass matrix term
184  if (flag == 2)
185  {
186  mass_matrix(local_eqn, local_unknown) +=
187  scaled_peclet_st * test(l) * psi(l2) * dS * W;
188  }
189 
190  // Assemble Jacobian term
191  jacobian(local_eqn, local_unknown) -=
192  // radial terms
193  (dS * dpsidx(l2, 0) *
194  (scaled_peclet * wind[0] * test(l) + dtestdx(l, 0)) +
195  // azimuthal terms
196  (sin_th * dpsidx(l2, 1) *
197  (r * scaled_peclet * wind[1] * test(l) + dtestdx(l, 1)))) *
198  W;
199  }
200  }
201  }
202  }
203  }
204 
205 
206  } // End of loop over integration points
207  }
208 
209 
210  //======================================================================
211  /// Self-test: Return 0 for OK
212  //======================================================================
213  // template <unsigned DIM>
215  {
216  bool passed = true;
217 
218  // Check lower-level stuff
219  if (FiniteElement::self_test() != 0)
220  {
221  passed = false;
222  }
223 
224  // Return verdict
225  if (passed)
226  {
227  return 0;
228  }
229  else
230  {
231  return 1;
232  }
233  }
234 
235 
236  //======================================================================
237  /// Output function:
238  ///
239  /// r,z,u,w_r,w_z
240  ///
241  /// nplot points in each coordinate direction
242  //======================================================================
244  const unsigned& nplot)
245  {
246  // Vector of local coordinates
247  Vector<double> s(2);
248 
249  // Tecplot header info
250  outfile << tecplot_zone_string(nplot);
251 
252  // Loop over plot points
253  unsigned num_plot_points = nplot_points(nplot);
254  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
255  {
256  // Get local coordinates of plot point
257  get_s_plot(iplot, nplot, s);
258 
259  // Get Eulerian coordinate of plot point
260  Vector<double> x(2);
261  interpolated_x(s, x);
262 
263  for (unsigned i = 0; i < 2; i++)
264  {
265  outfile << x[i] << " ";
266  }
267 
268  // Convert to Cartesians
269  outfile << x[0] * sin(x[1]) << " " << x[0] * cos(x[1]) << " ";
270  // Output concentration
271  outfile << this->interpolated_u_spherical_adv_diff(s) << " ";
272 
273  // Get the wind
274  Vector<double> wind(2);
275  // Dummy ipt argument needed... ?
276  unsigned ipt = 0;
277  get_wind_spherical_adv_diff(ipt, s, x, wind);
278  for (unsigned i = 0; i < 2; i++)
279  {
280  outfile << wind[i] << " ";
281  }
282  outfile << std::endl;
283  }
284 
285  // Write tecplot footer (e.g. FE connectivity lists)
286  write_tecplot_zone_footer(outfile, nplot);
287  }
288 
289 
290  //======================================================================
291  /// C-style output function:
292  ///
293  /// r,z,u
294  ///
295  /// nplot points in each coordinate direction
296  //======================================================================
297  // template <unsigned DIM>
299  const unsigned& nplot)
300  {
301  // Vector of local coordinates
302  Vector<double> s(2);
303 
304  // Tecplot header info
305  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
306 
307  // Loop over plot points
308  unsigned num_plot_points = nplot_points(nplot);
309  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
310  {
311  // Get local coordinates of plot point
312  get_s_plot(iplot, nplot, s);
313 
314  for (unsigned i = 0; i < 2; i++)
315  {
316  fprintf(file_pt, "%g ", interpolated_x(s, i));
317  }
318  // Write Cartesian coordinates
319  double r = interpolated_x(s, 0);
320  double theta = interpolated_x(s, 1);
321  fprintf(file_pt, "%g %g ", r * sin(theta), r * cos(theta));
322 
323  fprintf(file_pt, "%g \n", interpolated_u_spherical_adv_diff(s));
324  }
325 
326  // Write tecplot footer (e.g. FE connectivity lists)
327  write_tecplot_zone_footer(file_pt, nplot);
328  }
329 
330  //======================================================================
331  /// Output exact solution
332  ///
333  /// Solution is provided via function pointer.
334  /// Plot at a given number of plot points.
335  ///
336  /// r,z,u_exact
337  //======================================================================
338  // template <unsigned DIM>
340  std::ostream& outfile,
341  const unsigned& nplot,
343  {
344  // Vector of local coordinates
345  Vector<double> s(2);
346 
347  // Vector for coordintes
348  Vector<double> x(2);
349 
350  // Tecplot header info
351  outfile << tecplot_zone_string(nplot);
352 
353  // Exact solution Vector (here a scalar)
354  Vector<double> exact_soln(1);
355 
356  // Loop over plot points
357  unsigned num_plot_points = nplot_points(nplot);
358  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
359  {
360  // Get local coordinates of plot point
361  get_s_plot(iplot, nplot, s);
362 
363  // Get x position as Vector
364  interpolated_x(s, x);
365 
366  // Get exact solution at this point
367  (*exact_soln_pt)(x, exact_soln);
368 
369  // Output x,y,...,u_exact
370  for (unsigned i = 0; i < 2; i++)
371  {
372  outfile << x[i] << " ";
373  }
374  outfile << exact_soln[0] << std::endl;
375  }
376 
377  // Write tecplot footer (e.g. FE connectivity lists)
378  write_tecplot_zone_footer(outfile, nplot);
379  }
380 
381 
382  //======================================================================
383  /// Validate against exact solution
384  ///
385  /// Solution is provided via function pointer.
386  /// Plot error at a given number of plot points.
387  ///
388  //======================================================================
389  // template <unsigned DIM>
391  std::ostream& outfile,
393  double& error,
394  double& norm)
395  {
396  // Initialise
397  error = 0.0;
398  norm = 0.0;
399 
400  // Vector of local coordinates
401  Vector<double> s(2);
402 
403  // Vector for coordintes
404  Vector<double> x(2);
405 
406  // Find out how many nodes there are in the element
407  unsigned n_node = nnode();
408 
409  Shape psi(n_node);
410 
411  // Set the value of n_intpt
412  unsigned n_intpt = integral_pt()->nweight();
413 
414  // Tecplot header info
415  outfile << "ZONE" << std::endl;
416 
417  // Exact solution Vector (here a scalar)
418  Vector<double> exact_soln(1);
419 
420  // Loop over the integration points
421  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
422  {
423  // Assign values of s
424  for (unsigned i = 0; i < 2; i++)
425  {
426  s[i] = integral_pt()->knot(ipt, i);
427  }
428 
429  // Get the integral weight
430  double w = integral_pt()->weight(ipt);
431 
432  // Get jacobian of mapping
433  double J = J_eulerian(s);
434 
435  // Premultiply the weights and the Jacobian
436  double W = w * J;
437 
438  // Get x position as Vector
439  interpolated_x(s, x);
440 
441  // Get FE function value
442  double u_fe = interpolated_u_spherical_adv_diff(s);
443 
444  // Get exact solution at this point
445  (*exact_soln_pt)(x, exact_soln);
446 
447  // Output x,y,...,error
448  for (unsigned i = 0; i < 2; i++)
449  {
450  outfile << x[i] << " ";
451  }
452  outfile << exact_soln[0] << " " << exact_soln[0] - u_fe << std::endl;
453 
454  // Add to error and norm
455  norm += exact_soln[0] * exact_soln[0] * W;
456  error += (exact_soln[0] - u_fe) * (exact_soln[0] - u_fe) * W;
457  }
458  }
459 
460 
461  //======================================================================
462  // Set the data for the number of Variables at each node
463  //======================================================================
464  template<unsigned NNODE_1D>
466  1;
467 
468  //====================================================================
469  // Force build of templates
470  //====================================================================
471 
475 
476 
477 } // 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
TimeStepper *& time_stepper_pt()
Return the pointer to the timestepper.
Definition: nodes.h:238
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
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
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 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 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
virtual void get_source_spherical_adv_diff(const unsigned &ipt, const Vector< double > &x, double &source) const
Get source term at (Eulerian) position x. This function is virtual to allow overloading in multi-phys...
double interpolated_u_spherical_adv_diff(const Vector< double > &s) const
Return FE representation of function value u(s) at local coordinate s.
virtual double dshape_and_dtest_eulerian_at_knot_spherical_adv_diff(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Shape/test functions and derivs w.r.t. to global coords at integration point ipt; return Jacobian of ...
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact soln: r,z,u_exact at nplot^2 plot points.
void output(std::ostream &outfile)
Output with default number of plot points.
double du_dt_spherical_adv_diff(const unsigned &n) const
du/dt at local node n. Uses suitably interpolated value for hanging nodes.
const double & pe_st() const
Peclet number multiplied by Strouhal number.
void compute_error(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, double &error, double &norm)
Get error against and norm of exact solution.
virtual unsigned u_index_spherical_adv_diff() const
Return the index at which the unknown value is stored. The default value, 0, is appropriate for singl...
static double Default_peclet_number
Static default value for the Peclet number.
virtual void fill_in_generic_residual_contribution_spherical_adv_diff(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Add the element's contribution to its residual vector only (if flag=and/or element Jacobian matrix.
virtual void get_wind_spherical_adv_diff(const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &wind) const
Get wind at (Eulerian) position x and/or local coordinate s. This function is virtual to allow overlo...
virtual double weight(const unsigned &i, const unsigned &j) const
Access function for j-th weight for the i-th derivative.
Definition: timesteppers.h:594
//////////////////////////////////////////////////////////////////// ////////////////////////////////...