spherical_navier_stokes_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-2023 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 NS elements
27 
29 
30 
31 namespace oomph
32 {
33  /// Navier--Stokes equations static data
34  Vector<double> SphericalNavierStokesEquations::Gamma(3, 1.0);
35 
36  //=================================================================
37  /// "Magic" negative number that indicates that the pressure is
38  /// not stored at a node. This cannot be -1 because that represents
39  /// the positional hanging scheme in the hanging_pt object of nodes
40  //=================================================================
42 
43  /// Navier--Stokes equations static data
45 
46  /// Navier--Stokes equations static data
48 
49  /// Navier-Stokes equations default gravity vector
51 
52 
53  //================================================================
54  /// Compute the diagonal of the velocity/pressure mass matrices.
55  /// If which one=0, both are computed, otherwise only the pressure
56  /// (which_one=1) or the velocity mass matrix (which_one=2 -- the
57  /// LSC version of the preconditioner only needs that one)
58  /// NOTE: pressure versions isn't implemented yet because this
59  /// element has never been tried with Fp preconditoner.
60  //================================================================
63  Vector<double>& press_mass_diag,
64  Vector<double>& veloc_mass_diag,
65  const unsigned& which_one)
66  {
67 #ifdef PARANOID
68  if ((which_one == 0) || (which_one == 1))
69  {
70  throw OomphLibError("Computation of diagonal of pressure mass matrix is "
71  "not impmented yet.\n",
72  OOMPH_CURRENT_FUNCTION,
73  OOMPH_EXCEPTION_LOCATION);
74  }
75 #endif
76 
77  // Resize and initialise
78  veloc_mass_diag.assign(ndof(), 0.0);
79 
80  // find out how many nodes there are
81  const unsigned n_node = nnode();
82 
83  // find number of coordinates
84  const unsigned n_value = 3;
85 
86  // find the indices at which the local velocities are stored
87  Vector<unsigned> u_nodal_index(n_value);
88  for (unsigned i = 0; i < n_value; i++)
89  {
90  u_nodal_index[i] = this->u_index_spherical_nst(i);
91  }
92 
93  // Set up memory for test functions
94  Shape test(n_node);
95 
96  // Number of integration points
97  unsigned n_intpt = integral_pt()->nweight();
98 
99  // Integer to store the local equations no
100  int local_eqn = 0;
101 
102  // Loop over the integration points
103  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
104  {
105  // Get the integral weight
106  double w = integral_pt()->weight(ipt);
107 
108  // Get determinant of Jacobian of the mapping
109  double J = J_eulerian_at_knot(ipt);
110 
111  // Premultiply weights and Jacobian
112  double W = w * J;
113 
114  // Get the velocity test functions - there is no explicit
115  // function to give the test function so use shape
116  shape_at_knot(ipt, test);
117 
118  // Need to get the position to sort out the jacobian properly
119  double r = 0.0;
120  double theta = 0.0;
121  for (unsigned l = 0; l < n_node; l++)
122  {
123  r += this->nodal_position(l, 0) * test(l);
124  theta += this->nodal_position(l, 1) * test(l);
125  }
126 
127  // Multiply by the geometric factor
128  W *= r * r * sin(theta);
129 
130  // Loop over the veclocity test functions
131  for (unsigned l = 0; l < n_node; l++)
132  {
133  // Loop over the velocity components
134  for (unsigned i = 0; i < n_value; i++)
135  {
136  local_eqn = nodal_local_eqn(l, u_nodal_index[i]);
137 
138  // If not a boundary condition
139  if (local_eqn >= 0)
140  {
141  // Add the contribution
142  veloc_mass_diag[local_eqn] += test[l] * test[l] * W;
143  } // End of if not boundary condition statement
144  } // End of loop over dimension
145  } // End of loop over test functions
146  }
147  }
148 
149 
150  //======================================================================
151  /// Validate against exact velocity solution at given time.
152  /// Solution is provided via function pointer.
153  /// Plot at a given number of plot points and compute L2 error
154  /// and L2 norm of velocity solution over element.
155  //=======================================================================
157  std::ostream& outfile,
159  const double& time,
160  double& error,
161  double& norm)
162  {
163  error = 0.0;
164  norm = 0.0;
165 
166  // Vector of local coordinates
167  Vector<double> s(2);
168 
169  // Vector for coordintes
170  Vector<double> x(2);
171 
172  // Set the value of n_intpt
173  unsigned n_intpt = integral_pt()->nweight();
174 
175  outfile << "ZONE" << std::endl;
176 
177  // Exact solution Vector (u,v,[w],p)
178  Vector<double> exact_soln(4);
179 
180  // Loop over the integration points
181  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
182  {
183  // Assign values of s
184  for (unsigned i = 0; i < 2; i++)
185  {
186  s[i] = integral_pt()->knot(ipt, i);
187  }
188 
189  // Get the integral weight
190  double w = integral_pt()->weight(ipt);
191 
192  // Get jacobian of mapping
193  double J = J_eulerian(s);
194 
195  // Premultiply the weights and the Jacobian
196  double W = w * J;
197 
198  // Get x position as Vector
199  interpolated_x(s, x);
200 
201  // Get exact solution at this point
202  (*exact_soln_pt)(time, x, exact_soln);
203 
204  // Velocity error
205  for (unsigned i = 0; i < 3; i++)
206  {
207  norm += exact_soln[i] * exact_soln[i] * W;
208  error += (exact_soln[i] - interpolated_u_spherical_nst(s, i)) *
209  (exact_soln[i] - interpolated_u_spherical_nst(s, i)) * W;
210  }
211 
212  // Output x,y,...,u_exact
213  for (unsigned i = 0; i < 2; i++)
214  {
215  outfile << x[i] << " ";
216  }
217 
218  // Output x,y,[z],u_error,v_error,[w_error]
219  for (unsigned i = 0; i < 3; i++)
220  {
221  outfile << exact_soln[i] - interpolated_u_spherical_nst(s, i) << " ";
222  }
223 
224  outfile << std::endl;
225  }
226  }
227 
228 
229  //======================================================================
230  /// Validate against exact velocity solution
231  /// Solution is provided via function pointer.
232  /// Plot at a given number of plot points and compute L2 error
233  /// and L2 norm of velocity solution over element.
234  //=======================================================================
236  std::ostream& outfile,
238  double& error,
239  double& norm)
240  {
241  error = 0.0;
242  norm = 0.0;
243 
244  // Vector of local coordinates
245  Vector<double> s(2);
246 
247  // Vector for coordintes
248  Vector<double> x(2);
249 
250  // Set the value of n_intpt
251  unsigned n_intpt = integral_pt()->nweight();
252 
253 
254  outfile << "ZONE" << std::endl;
255 
256  // Exact solution Vector (u,v,w,p)
257  Vector<double> exact_soln(4);
258 
259  // Loop over the integration points
260  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
261  {
262  // Assign values of s
263  for (unsigned i = 0; i < 2; i++)
264  {
265  s[i] = integral_pt()->knot(ipt, i);
266  }
267 
268  // Get the integral weight
269  double w = integral_pt()->weight(ipt);
270 
271  // Get jacobian of mapping
272  double J = J_eulerian(s);
273 
274  // Premultiply the weights and the Jacobian
275  double W = w * J;
276 
277  // Get x position as Vector
278  interpolated_x(s, x);
279 
280  // Get the exact solution at this point
281  (*exact_soln_pt)(x, exact_soln);
282 
283  // Velocity error
284  for (unsigned i = 0; i < 3; i++)
285  {
286  norm += exact_soln[i] * exact_soln[i] * W;
287  error += (exact_soln[i] - interpolated_u_spherical_nst(s, i)) *
288  (exact_soln[i] - interpolated_u_spherical_nst(s, i)) * W;
289  }
290 
291  // Output x,y,...,u_exact
292  for (unsigned i = 0; i < 2; i++)
293  {
294  outfile << x[i] << " ";
295  }
296 
297  // Output x,y,[z],u_error,v_error,[w_error]
298  for (unsigned i = 0; i < 3; i++)
299  {
300  outfile << exact_soln[i] - interpolated_u_spherical_nst(s, i) << " ";
301  }
302  outfile << std::endl;
303  }
304  }
305 
306  //======================================================================
307  /// Validate against exact velocity solution
308  /// Solution is provided via direct coding.
309  /// Plot at a given number of plot points and compute energy error
310  /// and energy norm of velocity solution over element.
311  //=======================================================================
313  std::ostream& outfile,
316  FiniteElement::SteadyExactSolutionFctPt exact_soln_dtheta_pt,
317  double& u_error,
318  double& u_norm,
319  double& p_error,
320  double& p_norm)
321  {
322  u_error = 0.0;
323  p_error = 0.0;
324  u_norm = 0.0;
325  p_norm = 0.0;
326 
327  // Vector of local coordinates
328  Vector<double> s(2);
329 
330  // Vector for coordinates
331  Vector<double> x(2);
332 
333  // Set the value of n_intpt
334  unsigned n_intpt = integral_pt()->nweight();
335 
336 
337  /// Exact solution Vectors (u,v,[w],p) -
338  /// - need two extras to deal with the derivatives in the energy norm
339  Vector<double> exact_soln(4);
340  Vector<double> exact_soln_dr(4);
341  Vector<double> exact_soln_dth(4);
342 
343 
344  // Loop over the integration points
345  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
346  {
347  // Assign values of s
348  for (unsigned i = 0; i < 2; i++)
349  {
350  s[i] = integral_pt()->knot(ipt, i);
351  }
352 
353  // Get the integral weight
354  double w = integral_pt()->weight(ipt);
355 
356  // Get jacobian of mapping
357  double J = J_eulerian(s);
358 
359  // Premultiply the weights and the Jacobian
360  double W = w * J;
361 
362  // Get x position as Vector
363  interpolated_x(s, x);
364 
365  // Get exact solution at the point x using the 'actual' functions
366  (*exact_soln_pt)(x, exact_soln);
367  (*exact_soln_dr_pt)(x, exact_soln_dr);
368  (*exact_soln_dtheta_pt)(x, exact_soln_dth);
369 
370  // exact_soln = actual(x);
371  // exact_soln_dr = actual_dr(x);
372  // exact_soln_dth = actual_dth(x);
373 
374 
375  double r = x[0];
376  double theta = x[1];
377 
378  W *= r * r * sin(theta);
379 
380 
381  // Velocity error in energy norm
382  // Owing to the additional terms arising from the phi-components, we can't
383  // loop over the velocities. Calculate each section separately instead.
384 
385  // Norm calculation involves the double contraction of two tensors, so we
386  // take the summation of nine separate components for clarity. Fortunately
387  // most of the components zero out in the axisymmetric case.
388  //---------------------------------------------
389  // Entry (1,1) in dyad matrix
390  // No contribution
391  // Entry (1,2) in dyad matrix
392  // No contribution
393  // Entry (1,3) in dyad matrix
394  u_norm += (1 / (r * r)) * exact_soln[2] * exact_soln[2] * W;
395  //----------------------------------------------
396  // Entry (2,1) in dyad matrix
397  // No contribution
398  // Entry (2,2) in dyad matrix
399  // No contribution
400  // Entry (2,3) in dyad matrix
401  u_norm += cot(theta) * cot(theta) * (1 / (r * r)) * exact_soln[2] *
402  exact_soln[2] * W;
403  //-----------------------------------------------
404  // Entry (3,1) in dyad matrix
405  u_norm += exact_soln_dr[2] * exact_soln_dr[2] * W;
406  // Entry (3,2) in dyad matrix
407  u_norm += (1 / (r * r)) * exact_soln_dth[2] * exact_soln_dth[2] * W;
408  // Entry (3,3) in dyad matrix
409  // No contribution
410  //-----------------------------------------------
411  // End of norm calculation
412 
413 
414  // Error calculation involves the summation of 9 squared components,
415  // stated separately here for clarity
416  //---------------------------------------------
417  // Entry (1,1) in dyad matrix
418  u_error += interpolated_dudx_spherical_nst(s, 0, 0) *
420  // Entry (1,2) in dyad matrix
421  u_error += (1 / (r * r)) *
426  W;
427  // Entry (1,3) in dyad matrix
428  u_error += (1 / (r * r)) *
429  (interpolated_u_spherical_nst(s, 2) - exact_soln[2]) *
430  (interpolated_u_spherical_nst(s, 2) - exact_soln[2]) * W;
431  //---------------------------------------------
432  // Entry (2,1) in dyad matrix
433  u_error += interpolated_dudx_spherical_nst(s, 1, 0) *
435  // Entry (2,2) in dyad matrix
436  u_error += (1 / (r * r)) *
441  W;
442  // Entry (2,3) in dyad matrix
443  u_error += (1 / (r * r)) * cot(theta) * cot(theta) *
444  (interpolated_u_spherical_nst(s, 2) - exact_soln[2]) *
445  (interpolated_u_spherical_nst(s, 2) - exact_soln[2]) * W;
446  //---------------------------------------------
447  // Entry (3,1) in dyad matrix
448  u_error += (exact_soln_dr[2] - interpolated_dudx_spherical_nst(s, 2, 0)) *
449  (exact_soln_dr[2] - interpolated_dudx_spherical_nst(s, 2, 0)) *
450  W;
451  // Entry (3,2) in dyad matrix
452  u_error +=
453  (1 / (r * r)) *
454  (exact_soln_dth[2] - interpolated_dudx_spherical_nst(s, 2, 1)) *
455  (exact_soln_dth[2] - interpolated_dudx_spherical_nst(s, 2, 1)) * W;
456  // Entry (3,3) in dyad matrix
457  u_error += (1 / (r * r)) *
459  cot(theta) * interpolated_u_spherical_nst(s, 1)) *
461  cot(theta) * interpolated_u_spherical_nst(s, 1)) *
462  W;
463  //---------------------------------------------
464  // End of velocity error calculation
465 
466  // Pressure error in L2 norm - energy norm not required here as the
467  // pressure only appears undifferentiated in the NS SPC weak form.
468 
469  p_norm += exact_soln[3] * exact_soln[3] * W;
470 
471  p_error += (exact_soln[3] - interpolated_p_spherical_nst(s)) *
472  (exact_soln[3] - interpolated_p_spherical_nst(s)) * W;
473 
474 
475  // Output r and theta coordinates
476  for (unsigned i = 0; i < 2; i++)
477  {
478  outfile << x[i] << " ";
479  }
480 
481  // Output the velocity error details in the energy norm
482  outfile << u_error << " " << u_norm << " ";
483  // Output the pressure error details in the L2 norm
484  outfile << p_error << " " << p_norm << " ";
485  outfile << std::endl;
486 
487  } // End loop over integration points
488 
489  } // End of function statement
490 
491 
492  //======================================================================
493  // Output the shear stress at the outer wall of a rotating chamber.
494  //======================================================================
496  std::ostream& outfile)
497  {
498  // Vector of local coordinates
499  Vector<double> s(2);
500 
501  // Vector for coordinates
502  Vector<double> x(2);
503 
504  // Number of plot points
505  unsigned Npts = 3;
506 
507  // Declaration of the shear stress variable
508  double shear = 0.0;
509 
510 
511  // Assign values of s
512  for (unsigned i = 0; i < Npts; i++)
513  {
514  s[0] = 1.0;
515  s[1] = -1 + 2.0 * i / (Npts - 1);
516 
517 
518  // Get position as global coordinate
519  interpolated_x(s, x);
520 
521  // Output r and theta coordinates
522  outfile << x[0] << " ";
523  outfile << x[1] << " ";
524 
525  // Output the shear stress at the current coordinate
526  double r = x[0];
527  shear = interpolated_dudx_spherical_nst(s, 2, 0) -
528  (1 / r) * interpolated_u_spherical_nst(s, 2);
529  outfile << shear << " ";
530 
531  outfile << std::endl;
532  }
533  }
534 
535 
536  //======================================================================
537  // Output given velocity values and shear stress along a specific
538  // section of a rotating chamber.
539  //======================================================================
541  {
542  // Vector of local coordinates
543  Vector<double> s(2);
544 
545  // Vector for coordinates
546  Vector<double> x(2);
547 
548  // Number of plot points
549  unsigned Npts = 3;
550 
551 
552  // Assign values of s
553  for (unsigned i = 0; i < Npts; i++)
554  {
555  s[1] = 1.0;
556  s[0] = -1 + 2.0 * i / (Npts - 1);
557 
558 
559  // Get position as global coordinate
560  interpolated_x(s, x);
561  double r = x[0];
562  double theta = x[1];
563 
564  // Output r and theta coordinates
565  outfile << r << " ";
566  outfile << theta << " ";
567 
568  // Output the velocity at the current coordinate
569  outfile << interpolated_u_spherical_nst(s, 0) << " ";
570  outfile << interpolated_u_spherical_nst(s, 1) << " ";
571  outfile << interpolated_u_spherical_nst(s, 2) << " ";
572  outfile << interpolated_p_spherical_nst(s) << " ";
573 
574  // Output shear stress along the given radius
575  double shear = interpolated_dudx_spherical_nst(s, 2, 0) -
576  (1 / r) * interpolated_u_spherical_nst(s, 2);
577  outfile << shear << " ";
578 
579  // Coordinate setup for norm
580  double J = J_eulerian(s);
581  if (i != 1)
582  {
583  double w = integral_pt()->weight(i);
584  double W = w * J;
585  W *= r * r * sin(theta);
586 
587  // Output the velocity integral at current point
588  double u_int = 0.0;
589  for (unsigned j = 0; j < 3; j++)
590  {
591  u_int += interpolated_u_spherical_nst(s, j) * W;
592  }
593  outfile << u_int << " ";
594 
595  // Output the pressure integral at current point
596  double p_int = interpolated_p_spherical_nst(s) * W;
597  outfile << p_int;
598  outfile << std::endl;
599  }
600 
601  else
602  {
603  double w = integral_pt()->weight(5);
604  double W = w * J;
605  W *= r * r * sin(theta);
606 
607  // Output the velocity integral at current point
608  double u_int = 0.0;
609  for (unsigned j = 0; j < 3; j++)
610  {
611  u_int += interpolated_u_spherical_nst(s, j) * W;
612  }
613  outfile << u_int << " ";
614 
615  // Output the pressure integral at current point
616  double p_int = interpolated_p_spherical_nst(s) * W;
617  outfile << p_int;
618  outfile << std::endl;
619  }
620  }
621  }
622 
623 
624  //======================================================================
625  /// Output "exact" solution
626  /// Solution is provided via function pointer.
627  /// Plot at a given number of plot points.
628  /// Function prints as many components as are returned in solution Vector.
629  //=======================================================================
631  std::ostream& outfile,
632  const unsigned& nplot,
634  {
635  // Vector of local coordinates
636  Vector<double> s(2);
637 
638  // Vector for coordintes
639  Vector<double> x(2);
640 
641  // Tecplot header info
642  outfile << tecplot_zone_string(nplot);
643 
644  // Exact solution Vector
645  Vector<double> exact_soln;
646 
647 
648  // Loop over plot points
649  unsigned num_plot_points = nplot_points(nplot);
650  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
651  {
652  // Get local coordinates of plot point
653  get_s_plot(iplot, nplot, s);
654 
655  // Get x position as Vector
656  interpolated_x(s, x);
657 
658  // Get exact solution at this point
659  (*exact_soln_pt)(x, exact_soln);
660 
661  // Output x,y,...
662  for (unsigned i = 0; i < 3; i++)
663  {
664  outfile << x[i] << " ";
665  }
666 
667  // Output "exact solution"
668  for (unsigned i = 0; i < exact_soln.size(); i++)
669  {
670  outfile << exact_soln[i] << " ";
671  }
672 
673  outfile << std::endl;
674  }
675 
676  // Write tecplot footer (e.g. FE connectivity lists)
677  write_tecplot_zone_footer(outfile, nplot);
678  }
679 
680  //======================================================================
681  /// Output "exact" solution at a given time
682  /// Solution is provided via function pointer.
683  /// Plot at a given number of plot points.
684  /// Function prints as many components as are returned in solution Vector.
685  //=======================================================================
687  std::ostream& outfile,
688  const unsigned& nplot,
689  const double& time,
691  {
692  // Vector of local coordinates
693  Vector<double> s(2);
694 
695  // Vector for coordinates
696  Vector<double> x(2);
697 
698  // Tecplot header info
699  outfile << tecplot_zone_string(nplot);
700 
701  // Exact solution Vector
702  Vector<double> exact_soln;
703 
704  // Loop over plot points
705  unsigned num_plot_points = nplot_points(nplot);
706  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
707  {
708  // Get local coordinates of plot point
709  get_s_plot(iplot, nplot, s);
710 
711  // Get x position as Vector
712  interpolated_x(s, x);
713 
714  // Get exact solution at this point
715  (*exact_soln_pt)(time, x, exact_soln);
716 
717  // Output x,y,...
718  for (unsigned i = 0; i < 3; i++)
719  {
720  outfile << x[i] << " ";
721  }
722 
723  // Output "exact solution"
724  for (unsigned i = 0; i < exact_soln.size(); i++)
725  {
726  outfile << exact_soln[i] << " ";
727  }
728 
729  outfile << std::endl;
730  }
731 
732  // Write tecplot footer (e.g. FE connectivity lists)
733  write_tecplot_zone_footer(outfile, nplot);
734  }
735 
736  //==============================================================
737  /// Output function: Velocities only
738  /// x,y,[z],u,v,[w]
739  /// in tecplot format at specified previous timestep (t=0: present;
740  /// t>0: previous timestep). Specified number of plot points in each
741  /// coordinate direction.
742  //==============================================================
744  const unsigned& nplot,
745  const unsigned& t)
746  {
747  // Find number of nodes
748  unsigned n_node = nnode();
749 
750  // Local shape function
751  Shape psi(n_node);
752 
753  // Vectors of local coordinates and coords and velocities
754  Vector<double> s(2);
756  Vector<double> interpolated_u(3);
757 
758  // Tecplot header info
759  outfile << tecplot_zone_string(nplot);
760 
761  // Loop over plot points
762  unsigned num_plot_points = nplot_points(nplot);
763  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
764  {
765  // Get local coordinates of plot point
766  get_s_plot(iplot, nplot, s);
767 
768  // Get shape functions
769  shape(s, psi);
770 
771  // Loop over directions
772  for (unsigned i = 0; i < 3; i++)
773  {
774  interpolated_x[i] = 0.0;
775  interpolated_u[i] = 0.0;
776  // Get the index at which velocity i is stored
777  unsigned u_nodal_index = u_index_spherical_nst(i);
778  // Loop over the local nodes and sum
779  for (unsigned l = 0; l < n_node; l++)
780  {
781  interpolated_u[i] += nodal_value(t, l, u_nodal_index) * psi[l];
782  interpolated_x[i] += nodal_position(t, l, i) * psi[l];
783  }
784  }
785 
786  // Coordinates
787  for (unsigned i = 0; i < 2; i++)
788  {
789  outfile << interpolated_x[i] << " ";
790  }
791 
792  // Velocities
793  for (unsigned i = 0; i < 3; i++)
794  {
795  outfile << interpolated_u[i] << " ";
796  }
797 
798  outfile << std::endl;
799  }
800 
801  // Write tecplot footer (e.g. FE connectivity lists)
802  write_tecplot_zone_footer(outfile, nplot);
803  }
804 
805  //==============================================================
806  /// Output function:
807  /// x,y,[z],u,v,[w],p
808  /// in tecplot format. Specified number of plot points in each
809  /// coordinate direction.
810  //==============================================================
811  void SphericalNavierStokesEquations::output(std::ostream& outfile,
812  const unsigned& nplot)
813  {
814  // Vector of local coordinates
815  Vector<double> s(2);
816 
817  // Tecplot header info
818  outfile << tecplot_zone_string(nplot);
819 
820  // Loop over plot points
821  unsigned num_plot_points = nplot_points(nplot);
822  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
823  {
824  // Get local coordinates of plot point
825  get_s_plot(iplot, nplot, s);
826 
827  // Coordinates in Cartesian Coordinates
828  // for(unsigned i=0;i<DIM;i++)
829  // {
830  // outfile << interpolated_x(s,i) << " ";
831  // }
832 
833  // Output the coordinates in Cartesian coordintes
834  double r = interpolated_x(s, 0);
835  double theta = interpolated_x(s, 1);
836  // The actual Cartesian values
837  outfile << r * sin(theta) << " " << r * cos(theta) << " ";
838 
839  // Now the x and y velocities
840  double u_r = interpolated_u_spherical_nst(s, 0);
841  double u_theta = interpolated_u_spherical_nst(s, 1);
842 
843  outfile << u_r * sin(theta) + u_theta * cos(theta) << " "
844  << u_r * cos(theta) - u_theta * sin(theta) << " ";
845 
846  // Now the polar coordinates
847  // outfile << r << " " << theta << " ";
848 
849  // Velocities r, theta, phi
850  for (unsigned i = 0; i < 3; i++)
851  {
852  outfile << interpolated_u_spherical_nst(s, i) << " ";
853  }
854 
855  // Pressure
856  outfile << interpolated_p_spherical_nst(s) << " ";
857 
858  // Dissipation
859  outfile << dissipation(s) << " ";
860 
861  outfile << std::endl;
862  }
863 
864  // Write tecplot footer (e.g. FE connectivity lists)
865  write_tecplot_zone_footer(outfile, nplot);
866  }
867 
868 
869  //==============================================================
870  /// C-style output function:
871  /// x,y,[z],u,v,[w],p
872  /// in tecplot format. Specified number of plot points in each
873  /// coordinate direction.
874  //==============================================================
876  const unsigned& nplot)
877  {
878  // Vector of local coordinates
879  Vector<double> s(2);
880 
881  // Tecplot header info
882  fprintf(file_pt, "%s", tecplot_zone_string(nplot).c_str());
883 
884  // Loop over plot points
885  unsigned num_plot_points = nplot_points(nplot);
886  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
887  {
888  // Get local coordinates of plot point
889  get_s_plot(iplot, nplot, s);
890 
891  // Coordinates
892  for (unsigned i = 0; i < 2; i++)
893  {
894  fprintf(file_pt, "%g ", interpolated_x(s, i));
895  }
896 
897  // Velocities
898  for (unsigned i = 0; i < 3; i++)
899  {
900  fprintf(file_pt, "%g ", interpolated_u_spherical_nst(s, i));
901  }
902 
903  // Pressure
904  fprintf(file_pt, "%g \n", interpolated_p_spherical_nst(s));
905  }
906  fprintf(file_pt, "\n");
907 
908  // Write tecplot footer (e.g. FE connectivity lists)
909  write_tecplot_zone_footer(file_pt, nplot);
910  }
911 
912 
913  //==============================================================
914  /// Full output function:
915  /// x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation
916  /// in tecplot format. Specified number of plot points in each
917  /// coordinate direction
918  //==============================================================
920  const unsigned& nplot)
921  {
922  throw OomphLibError(
923  "Probably Broken", OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
924 
925  // Vector of local coordinates
926  Vector<double> s(2);
927 
928  // Acceleration
929  Vector<double> dudt(3);
930 
931  // Mesh elocity
932  Vector<double> mesh_veloc(3, 0.0);
933 
934  // Tecplot header info
935  outfile << tecplot_zone_string(nplot);
936 
937  // Find out how many nodes there are
938  unsigned n_node = nnode();
939 
940  // Set up memory for the shape functions
941  Shape psif(n_node);
942  DShape dpsifdx(n_node, 2);
943 
944  // Loop over plot points
945  unsigned num_plot_points = nplot_points(nplot);
946  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
947  {
948  // Get local coordinates of plot point
949  get_s_plot(iplot, nplot, s);
950 
951  // Call the derivatives of the shape and test functions
952  dshape_eulerian(s, psif, dpsifdx);
953 
954  // Allocate storage
955  Vector<double> mesh_veloc(3);
956  Vector<double> dudt(3);
957  Vector<double> dudt_ALE(3);
958  DenseMatrix<double> interpolated_dudx(3, 2);
959 
960  // Initialise everything to zero
961  for (unsigned i = 0; i < 3; i++)
962  {
963  mesh_veloc[i] = 0.0;
964  dudt[i] = 0.0;
965  dudt_ALE[i] = 0.0;
966  for (unsigned j = 0; j < 2; j++)
967  {
968  interpolated_dudx(i, j) = 0.0;
969  }
970  }
971 
972  // Calculate velocities and derivatives
973 
974 
975  // Loop over directions
976  for (unsigned i = 0; i < 3; i++)
977  {
978  // Get the index at which velocity i is stored
979  unsigned u_nodal_index = u_index_spherical_nst(i);
980  // Loop over nodes
981  for (unsigned l = 0; l < n_node; l++)
982  {
983  dudt[i] += du_dt_spherical_nst(l, u_nodal_index) * psif(l);
984  mesh_veloc[i] += dnodal_position_dt(l, i) * psif(l);
985 
986  // Loop over derivative directions for velocity gradients
987  for (unsigned j = 0; j < 2; j++)
988  {
989  interpolated_dudx(i, j) +=
990  nodal_value(l, u_nodal_index) * dpsifdx(l, j);
991  }
992  }
993  }
994 
995 
996  // Get dudt in ALE form (incl mesh veloc)
997  for (unsigned i = 0; i < 3; i++)
998  {
999  dudt_ALE[i] = dudt[i];
1000  for (unsigned k = 0; k < 2; k++)
1001  {
1002  dudt_ALE[i] -= mesh_veloc[k] * interpolated_dudx(i, k);
1003  }
1004  }
1005 
1006 
1007  // Coordinates
1008  for (unsigned i = 0; i < 2; i++)
1009  {
1010  outfile << interpolated_x(s, i) << " ";
1011  }
1012 
1013  // Velocities
1014  for (unsigned i = 0; i < 3; i++)
1015  {
1016  outfile << interpolated_u_spherical_nst(s, i) << " ";
1017  }
1018 
1019  // Pressure
1020  outfile << interpolated_p_spherical_nst(s) << " ";
1021 
1022  // Accelerations
1023  for (unsigned i = 0; i < 3; i++)
1024  {
1025  outfile << dudt_ALE[i] << " ";
1026  }
1027 
1028  // Dissipation
1029  outfile << dissipation(s) << " ";
1030 
1031 
1032  outfile << std::endl;
1033  }
1034 
1035  // Write tecplot footer (e.g. FE connectivity lists)
1036  write_tecplot_zone_footer(outfile, nplot);
1037  }
1038 
1039 
1040  //==============================================================
1041  /// Output function for vorticity.
1042  /// x,y,[z],[omega_x,omega_y,[and/or omega_z]]
1043  /// in tecplot format. Specified number of plot points in each
1044  /// coordinate direction.
1045  //==============================================================
1047  const unsigned& nplot)
1048  {
1049  // Vector of local coordinates
1050  Vector<double> s(2);
1051 
1052  // Create vorticity vector of the required size
1053  Vector<double> vorticity(1);
1054 
1055  // Tecplot header info
1056  outfile << tecplot_zone_string(nplot);
1057 
1058  // Loop over plot points
1059  unsigned num_plot_points = nplot_points(nplot);
1060  for (unsigned iplot = 0; iplot < num_plot_points; iplot++)
1061  {
1062  // Get local coordinates of plot point
1063  get_s_plot(iplot, nplot, s);
1064 
1065  // Coordinates
1066  for (unsigned i = 0; i < 2; i++)
1067  {
1068  outfile << interpolated_x(s, i) << " ";
1069  }
1070 
1071  // Get vorticity
1072  get_vorticity(s, vorticity);
1073 
1074  outfile << vorticity[0] << std::endl;
1075  ;
1076  }
1077 
1078  // Write tecplot footer (e.g. FE connectivity lists)
1079  write_tecplot_zone_footer(outfile, nplot);
1080  }
1081 
1082 
1083  //==============================================================
1084  /// Return integral of dissipation over element
1085  //==============================================================
1087  {
1088  // Initialise
1089  double diss = 0.0;
1090 
1091  // Set the value of n_intpt
1092  unsigned n_intpt = integral_pt()->nweight();
1093 
1094  // Set the Vector to hold local coordinates
1095  Vector<double> s(2);
1096 
1097  // Loop over the integration points
1098  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1099  {
1100  // Assign values of s
1101  for (unsigned i = 0; i < 2; i++)
1102  {
1103  s[i] = integral_pt()->knot(ipt, i);
1104  }
1105 
1106  // Get the integral weight
1107  double w = integral_pt()->weight(ipt);
1108 
1109  // Get Jacobian of mapping
1110  double J = J_eulerian(s);
1111 
1112  // Get strain rate matrix
1113  DenseMatrix<double> strainrate(3, 3);
1114  strain_rate(s, strainrate);
1115 
1116  // Initialise
1117  double local_diss = 0.0;
1118  for (unsigned i = 0; i < 3; i++)
1119  {
1120  for (unsigned j = 0; j < 3; j++)
1121  {
1122  local_diss += 2.0 * strainrate(i, j) * strainrate(i, j);
1123  }
1124  }
1125 
1126  diss += local_diss * w * J;
1127  }
1128 
1129  return diss;
1130  }
1131 
1132  //==============================================================
1133  /// Compute traction (on the viscous scale) exerted onto
1134  /// the fluid at local coordinate s. N has to be outer unit normal
1135  /// to the fluid.
1136  //==============================================================
1138  const Vector<double>& N,
1139  Vector<double>& traction)
1140  {
1141  // Get velocity gradients
1142  DenseMatrix<double> strainrate(3, 3);
1143  strain_rate(s, strainrate);
1144 
1145  // Get pressure
1146  double press = interpolated_p_spherical_nst(s);
1147 
1148  // Loop over traction components
1149  for (unsigned i = 0; i < 3; i++)
1150  {
1151  // If we are in the r and theta direction
1152  // initialise the traction
1153  if (i < 2)
1154  {
1155  traction[i] = -press * N[i];
1156  }
1157  // Otherwise it's zero because the normal cannot have a component
1158  // in the phi direction
1159  else
1160  {
1161  traction[i] = 0.0;
1162  }
1163  // Loop over the possible traction directions
1164  for (unsigned j = 0; j < 2; j++)
1165  {
1166  traction[i] += 2.0 * strainrate(i, j) * N[j];
1167  }
1168  }
1169  }
1170 
1171  //==============================================================
1172  /// Return dissipation at local coordinate s
1173  //==============================================================
1175  const Vector<double>& s) const
1176  {
1177  // Get strain rate matrix
1178  DenseMatrix<double> strainrate(3, 3);
1179  strain_rate(s, strainrate);
1180 
1181  // Initialise
1182  double local_diss = 0.0;
1183  for (unsigned i = 0; i < 3; i++)
1184  {
1185  for (unsigned j = 0; j < 3; j++)
1186  {
1187  local_diss += 2.0 * strainrate(i, j) * strainrate(i, j);
1188  }
1189  }
1190 
1191  return local_diss;
1192  }
1193 
1194  //==============================================================
1195  /// Get strain-rate tensor: 1/2 (du_i/dx_j + du_j/dx_i)
1196  //==============================================================
1198  const Vector<double>& s, DenseMatrix<double>& strainrate) const
1199  {
1200  // Velocity components
1201  Vector<double> interpolated_u(3, 0.0);
1202  // Coordinates
1203  double interpolated_r = 0.0;
1204  double interpolated_theta = 0.0;
1205  // Velocity gradient matrix
1206  DenseMatrix<double> dudx(3, 2, 0.0);
1207 
1208 
1209  // Find out how many nodes there are in the element
1210  unsigned n_node = nnode();
1211 
1212  // Set up memory for the shape and test functions
1213  Shape psi(n_node);
1214  DShape dpsidx(n_node, 2);
1215 
1216  // Call the derivatives of the shape functions
1217  dshape_eulerian(s, psi, dpsidx);
1218 
1219  // Get the values of the positions
1220  for (unsigned l = 0; l < n_node; l++)
1221  {
1222  interpolated_r += this->nodal_position(l, 0) * psi(l);
1223  interpolated_theta += this->nodal_position(l, 1) * psi(l);
1224  }
1225 
1226 
1227  // Loop over velocity components
1228  for (unsigned i = 0; i < 3; i++)
1229  {
1230  // Get the index at which the i-th velocity is stored
1231  unsigned u_nodal_index = u_index_spherical_nst(i);
1232  // Loop over nodes
1233  for (unsigned l = 0; l < n_node; l++)
1234  {
1235  // Get the velocity value
1236  const double nodal_u = this->nodal_value(l, u_nodal_index);
1237  // Add in the velocities
1238  interpolated_u[i] += nodal_u * psi(l);
1239  // Loop over the derivative directions
1240  for (unsigned j = 0; j < 2; j++)
1241  {
1242  dudx(i, j) += nodal_u * dpsidx(l, j);
1243  }
1244  }
1245  }
1246 
1247  // Assign those strain rates without any negative powers of the radius
1248  strainrate(0, 0) = dudx(0, 0);
1249  strainrate(0, 1) = 0.0;
1250  strainrate(0, 2) = 0.0;
1251  ;
1252  strainrate(1, 1) = 0.0;
1253  strainrate(1, 2) = 0.0;
1254  strainrate(2, 2) = 0.0;
1255 
1256 
1257  // Add the negative powers of the radius, unless we are
1258  // at the origin
1259  if (std::fabs(interpolated_r) > 1.0e-15)
1260  {
1261  const double pi = MathematicalConstants::Pi;
1262  double inverse_r = 1.0 / interpolated_r;
1263 
1264  // Should we include the cot terms (default no)
1265  bool include_cot_terms = false;
1266  double cot_theta = 0.0;
1267  // If we in the legal range then include cot
1268  if ((std::fabs(interpolated_theta) > 1.0e-15) &&
1269  (std::fabs(pi - interpolated_theta) > 1.0e-15))
1270  {
1271  include_cot_terms = true;
1272  cot_theta = this->cot(interpolated_theta);
1273  }
1274 
1275  strainrate(0, 1) =
1276  0.5 * (dudx(1, 0) + (dudx(0, 1) - interpolated_u[1]) * inverse_r);
1277  strainrate(0, 2) = 0.5 * (dudx(2, 0) - interpolated_u[2] * inverse_r);
1278  strainrate(1, 1) = (dudx(1, 1) + interpolated_u[0]) * inverse_r;
1279 
1280  if (include_cot_terms)
1281  {
1282  strainrate(1, 2) =
1283  0.5 * (dudx(2, 1) - interpolated_u[2] * cot_theta) * inverse_r;
1284  strainrate(2, 2) =
1285  (interpolated_u[0] + interpolated_u[1] * cot_theta) * inverse_r;
1286  }
1287  }
1288 
1289  // Sort out the symmetries
1290  strainrate(1, 0) = strainrate(0, 1);
1291  strainrate(2, 0) = strainrate(0, 2);
1292  strainrate(2, 1) = strainrate(1, 2);
1293  }
1294 
1295 
1296  //==============================================================
1297  /// Compute 2D vorticity vector at local coordinate s (return in
1298  /// one and only component of vorticity vector
1299  //==============================================================
1301  const Vector<double>& s, Vector<double>& vorticity) const
1302  {
1303  throw OomphLibError("Not tested for spherical Navier-Stokes elements",
1304  OOMPH_CURRENT_FUNCTION,
1305  OOMPH_EXCEPTION_LOCATION);
1306 
1307 #ifdef PARANOID
1308  if (vorticity.size() != 1)
1309  {
1310  std::ostringstream error_message;
1311  error_message << "Computation of vorticity in 2D requires a 1D vector\n"
1312  << "which contains the only non-zero component of the\n"
1313  << "vorticity vector. You've passed a vector of size "
1314  << vorticity.size() << std::endl;
1315 
1316  throw OomphLibError(
1317  error_message.str(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION);
1318  }
1319 #endif
1320 
1321  // Specify spatial dimension
1322  unsigned DIM = 2;
1323 
1324  // Velocity gradient matrix
1325  DenseMatrix<double> dudx(DIM);
1326 
1327  // Find out how many nodes there are in the element
1328  unsigned n_node = nnode();
1329 
1330  // Set up memory for the shape and test functions
1331  Shape psi(n_node);
1332  DShape dpsidx(n_node, DIM);
1333 
1334  // Call the derivatives of the shape functions
1335  dshape_eulerian(s, psi, dpsidx);
1336 
1337  // Initialise to zero
1338  for (unsigned i = 0; i < DIM; i++)
1339  {
1340  for (unsigned j = 0; j < DIM; j++)
1341  {
1342  dudx(i, j) = 0.0;
1343  }
1344  }
1345 
1346  // Loop over veclocity components
1347  for (unsigned i = 0; i < DIM; i++)
1348  {
1349  // Get the index at which the i-th velocity is stored
1350  unsigned u_nodal_index = u_index_spherical_nst(i);
1351  // Loop over derivative directions
1352  for (unsigned j = 0; j < DIM; j++)
1353  {
1354  // Loop over nodes
1355  for (unsigned l = 0; l < n_node; l++)
1356  {
1357  dudx(i, j) += nodal_value(l, u_nodal_index) * dpsidx(l, j);
1358  }
1359  }
1360  }
1361 
1362  // Z-component of vorticity
1363  vorticity[0] = dudx(1, 0) - dudx(0, 1);
1364  }
1365 
1366 
1367  //==============================================================
1368  /// Get integral of kinetic energy over element:
1369  //==============================================================
1371  {
1372  // Initialise
1373  double kin_en = 0.0;
1374 
1375  // Set the value of n_intpt
1376  unsigned n_intpt = integral_pt()->nweight();
1377 
1378  // Set the Vector to hold local coordinates
1379  Vector<double> s(2);
1380 
1381  // Loop over the integration points
1382  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1383  {
1384  // Assign values of s
1385  for (unsigned i = 0; i < 2; i++)
1386  {
1387  s[i] = integral_pt()->knot(ipt, i);
1388  }
1389 
1390  // Get the integral weight
1391  double w = integral_pt()->weight(ipt);
1392 
1393  // Get Jacobian of mapping
1394  double J = J_eulerian(s);
1395 
1396  // Loop over directions
1397  double veloc_squared = 0.0;
1398  for (unsigned i = 0; i < 3; i++)
1399  {
1400  veloc_squared += interpolated_u_spherical_nst(s, i) *
1402  }
1403 
1404  kin_en += 0.5 * veloc_squared * w * J;
1405  }
1406 
1407  return kin_en;
1408  }
1409 
1410 
1411  //==========================================================================
1412  /// Get integral of time derivative of kinetic energy over element:
1413  //==========================================================================
1415  {
1416  throw OomphLibError("Not tested for spherical Navier-Stokes elements",
1417  OOMPH_CURRENT_FUNCTION,
1418  OOMPH_EXCEPTION_LOCATION);
1419 
1420 
1421  // Initialise
1422  double d_kin_en_dt = 0.0;
1423 
1424  // Set the value of n_intpt
1425  const unsigned n_intpt = integral_pt()->nweight();
1426 
1427  // Set the Vector to hold local coordinates
1428  Vector<double> s(2);
1429 
1430  // Get the number of nodes
1431  const unsigned n_node = this->nnode();
1432 
1433  // Storage for the shape function
1434  Shape psi(n_node);
1435  DShape dpsidx(n_node, 2);
1436 
1437  // Get the value at which the velocities are stored
1438  unsigned u_index[3];
1439  for (unsigned i = 0; i < 3; i++)
1440  {
1441  u_index[i] = this->u_index_spherical_nst(i);
1442  }
1443 
1444  // Loop over the integration points
1445  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1446  {
1447  // Get the jacobian of the mapping
1448  double J = dshape_eulerian_at_knot(ipt, psi, dpsidx);
1449 
1450  // Get the integral weight
1451  double w = integral_pt()->weight(ipt);
1452 
1453  // Now work out the velocity and the time derivative
1454  Vector<double> interpolated_u(3, 0.0), interpolated_dudt(3, 0.0);
1455 
1456  // Loop over the shape functions
1457  for (unsigned l = 0; l < n_node; l++)
1458  {
1459  // Loop over the velocity components
1460  for (unsigned i = 0; i < 3; i++)
1461  {
1462  interpolated_u[i] += nodal_value(l, u_index[i]) * psi(l);
1463  interpolated_dudt[i] += du_dt_spherical_nst(l, u_index[i]) * psi(l);
1464  }
1465  }
1466 
1467  // Get mesh velocity bit
1468  if (!ALE_is_disabled)
1469  {
1470  Vector<double> mesh_velocity(2, 0.0);
1471  DenseMatrix<double> interpolated_dudx(3, 2, 0.0);
1472 
1473  // Loop over nodes again
1474  for (unsigned l = 0; l < n_node; l++)
1475  {
1476  // Mesh can only move in the first two directions
1477  for (unsigned i = 0; i < 2; i++)
1478  {
1479  mesh_velocity[i] += this->dnodal_position_dt(l, i) * psi(l);
1480  }
1481 
1482  // Now du/dx bit
1483  for (unsigned i = 0; i < 3; i++)
1484  {
1485  for (unsigned j = 0; j < 2; j++)
1486  {
1487  interpolated_dudx(i, j) +=
1488  this->nodal_value(l, u_index[i]) * dpsidx(l, j);
1489  }
1490  }
1491  }
1492 
1493  // Subtract mesh velocity from du_dt
1494  for (unsigned i = 0; i < 3; i++)
1495  {
1496  for (unsigned k = 0; k < 2; k++)
1497  {
1498  interpolated_dudt[i] -= mesh_velocity[k] * interpolated_dudx(i, k);
1499  }
1500  }
1501  }
1502 
1503 
1504  // Loop over directions and add up u du/dt terms
1505  double sum = 0.0;
1506  for (unsigned i = 0; i < 3; i++)
1507  {
1508  sum += interpolated_u[i] * interpolated_dudt[i];
1509  }
1510 
1511  d_kin_en_dt += sum * w * J;
1512  }
1513 
1514  return d_kin_en_dt;
1515  }
1516 
1517 
1518  //==============================================================
1519  /// Return pressure integrated over the element
1520  //==============================================================
1522  {
1523  // Initialise
1524  double press_int = 0;
1525 
1526  // Set the value of n_intpt
1527  unsigned n_intpt = integral_pt()->nweight();
1528 
1529  // Set the Vector to hold local coordinates
1530  Vector<double> s(2);
1531 
1532  // Loop over the integration points
1533  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1534  {
1535  // Assign values of s
1536  for (unsigned i = 0; i < 2; i++)
1537  {
1538  s[i] = integral_pt()->knot(ipt, i);
1539  }
1540 
1541  // Get the integral weight
1542  double w = integral_pt()->weight(ipt);
1543 
1544  // Get Jacobian of mapping
1545  double J = J_eulerian(s);
1546 
1547  // Premultiply the weights and the Jacobian
1548  double W = w * J;
1549 
1550  // Get pressure
1551  double press = interpolated_p_spherical_nst(s);
1552 
1553  // Add
1554  press_int += press * W;
1555  }
1556 
1557  return press_int;
1558  }
1559 
1560  //==============================================================
1561  /// Compute the residuals for the Navier--Stokes
1562  /// equations; flag=1(or 0): do (or don't) compute the
1563  /// Jacobian as well.
1564  //==============================================================
1567  Vector<double>& residuals,
1568  DenseMatrix<double>& jacobian,
1569  DenseMatrix<double>& mass_matrix,
1570  unsigned flag)
1571  {
1572  // Return immediately if there are no dofs
1573  if (ndof() == 0) return;
1574 
1575  // Find out how many nodes there are
1576  const unsigned n_node = nnode();
1577 
1578  // Get continuous time from timestepper of first node
1579  double time = node_pt(0)->time_stepper_pt()->time_pt()->time();
1580 
1581  // Find out how many pressure dofs there are
1582  const unsigned n_pres = npres_spherical_nst();
1583 
1584  // Find the indices at which the local velocities are stored
1585  unsigned u_nodal_index[3];
1586  for (unsigned i = 0; i < 3; i++)
1587  {
1588  u_nodal_index[i] = u_index_spherical_nst(i);
1589  }
1590 
1591  // Set up memory for the shape and test functions
1592  Shape psif(n_node), testf(n_node);
1593  DShape dpsifdx(n_node, 2), dtestfdx(n_node, 2);
1594 
1595  // Set up memory for pressure shape and test functions
1596  Shape psip(n_pres), testp(n_pres);
1597 
1598  // Number of integration points
1599  const unsigned n_intpt = integral_pt()->nweight();
1600 
1601  // Set the Vector to hold local coordinates
1602  Vector<double> s(2);
1603 
1604  // Get Physical Variables from Element
1605  // Reynolds number must be multiplied by the density ratio
1606  const double scaled_re = re() * density_ratio();
1607  const double scaled_re_st = re_st() * density_ratio();
1608  const double scaled_re_inv_ro = re_invro() * density_ratio();
1609  // double scaled_re_inv_fr = re_invfr()*density_ratio();
1610  // double visc_ratio = viscosity_ratio();
1611  Vector<double> G = g();
1612 
1613  // Integers to store the local equations and unknowns
1614  int local_eqn = 0, local_unknown = 0;
1615 
1616  // Loop over the integration points
1617  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
1618  {
1619  // Assign values of s
1620  for (unsigned i = 0; i < 2; i++) s[i] = integral_pt()->knot(ipt, i);
1621  // Get the integral weight
1622  double w = integral_pt()->weight(ipt);
1623 
1624  // Call the derivatives of the shape and test functions
1626  ipt, psif, dpsifdx, testf, dtestfdx);
1627 
1628  // Call the pressure shape and test functions
1629  pshape_spherical_nst(s, psip, testp);
1630 
1631  // Premultiply the weights and the Jacobian
1632  const double W = w * J;
1633 
1634  // Calculate local values of the pressure and velocity components
1635  // Allocate
1636  double interpolated_p = 0.0;
1637  Vector<double> interpolated_u(3, 0.0);
1639  Vector<double> mesh_velocity(2, 0.0);
1640  Vector<double> dudt(3, 0.0);
1641  DenseMatrix<double> interpolated_dudx(3, 2, 0.0);
1642 
1643  // Calculate pressure
1644  for (unsigned l = 0; l < n_pres; l++)
1645  {
1646  interpolated_p += p_spherical_nst(l) * psip(l);
1647  }
1648 
1649 
1650  // Calculate velocities and derivatives:
1651 
1652  // Loop over nodes
1653  for (unsigned l = 0; l < n_node; l++)
1654  {
1655  // cache the shape functions
1656  double psi_ = psif(l);
1657  // Loop over the positions separately
1658  for (unsigned i = 0; i < 2; i++)
1659  {
1660  interpolated_x[i] += raw_nodal_position(l, i) * psi_;
1661  }
1662 
1663  // Loop over velocity directions (three of these)
1664  for (unsigned i = 0; i < 3; i++)
1665  {
1666  // Get the nodal value
1667  const double u_value = raw_nodal_value(l, u_nodal_index[i]);
1668  interpolated_u[i] += u_value * psi_;
1669  dudt[i] += du_dt_spherical_nst(l, i) * psi_;
1670 
1671  // Loop over derivative directions
1672  for (unsigned j = 0; j < 2; j++)
1673  {
1674  interpolated_dudx(i, j) += u_value * dpsifdx(l, j);
1675  }
1676  }
1677  }
1678 
1679  if (!ALE_is_disabled)
1680  {
1681  OomphLibWarning("ALE is not tested for spherical Navier Stokes",
1682  "SphericalNS::fill_in...",
1683  OOMPH_EXCEPTION_LOCATION);
1684  // Loop over nodes
1685  for (unsigned l = 0; l < n_node; l++)
1686  {
1687  // Loop over directions (only DIM (2) of these)
1688  for (unsigned i = 0; i < 2; i++)
1689  {
1690  mesh_velocity[i] += this->raw_dnodal_position_dt(l, i) * psif(l);
1691  }
1692  }
1693  }
1694 
1695  // Get the user-defined body force terms
1696  Vector<double> body_force(3);
1697  get_body_force_spherical_nst(time, ipt, s, interpolated_x, body_force);
1698 
1699  // Get the user-defined source function
1700  // double source = get_source_spherical_nst(time(),ipt,interpolated_x);
1701 
1702 
1703  // MOMENTUM EQUATIONS
1704  //------------------
1705 
1706  const double r = interpolated_x[0];
1707  const double r2 = r * r;
1708  const double sin_theta = sin(interpolated_x[1]);
1709  const double cos_theta = cos(interpolated_x[1]);
1710  const double cot_theta = cos_theta / sin_theta;
1711 
1712  const double u_r = interpolated_u[0];
1713  const double u_theta = interpolated_u[1];
1714  const double u_phi = interpolated_u[2];
1715 
1716  // Loop over the test functions
1717  for (unsigned l = 0; l < n_node; l++)
1718  {
1719  // Cannot loop over the velocity components,
1720  // so address each of them separately
1721 
1722  // FIRST: r-component momentum equation
1723 
1724  unsigned i = 0;
1725 
1726  // If it's not a boundary condition
1727  local_eqn = nodal_local_eqn(l, u_nodal_index[i]);
1728  if (local_eqn >= 0)
1729  {
1730  // Convective r-terms
1731  double conv = u_r * interpolated_dudx(0, 0) * r;
1732 
1733  // Convective theta-terms
1734  conv += u_theta * interpolated_dudx(0, 1);
1735 
1736  // Remaining convective terms
1737  conv -= (u_theta * u_theta + u_phi * u_phi);
1738 
1739  // Add the time-derivative term and the convective terms
1740  // to the sum
1741  double sum = (scaled_re_st * r2 * dudt[0] + r * scaled_re * conv) *
1742  sin_theta * testf[l];
1743 
1744  // Subtract the mesh velocity terms
1745  if (!ALE_is_disabled)
1746  {
1747  sum -= scaled_re_st *
1748  (mesh_velocity[0] * interpolated_dudx(0, 0) * r +
1749  mesh_velocity[1] * (interpolated_dudx(0, 1) - u_theta)) *
1750  r * sin_theta * testf[l];
1751  }
1752 
1753  // Add the Coriolis term
1754  sum -= 2.0 * scaled_re_inv_ro * sin_theta * u_phi * r2 * sin_theta *
1755  testf[l];
1756 
1757  // r-derivative test function component of stress tensor
1758  sum += (-interpolated_p + 2.0 * interpolated_dudx(0, 0)) * r2 *
1759  sin_theta * dtestfdx(l, 0);
1760 
1761  // theta-derivative test function component of stress tensor
1762  sum +=
1763  (r * interpolated_dudx(1, 0) - u_theta + interpolated_dudx(0, 1)) *
1764  sin_theta * dtestfdx(l, 1);
1765 
1766  // Undifferentiated test function component of stress tensor
1767  sum += 2.0 *
1768  ((-r * interpolated_p + interpolated_dudx(1, 1) + 2.0 * u_r) *
1769  sin_theta +
1770  u_theta * cos_theta) *
1771  testf(l);
1772 
1773  // Make the residuals negative for consistency with the
1774  // other Navier-Stokes equations
1775  residuals[local_eqn] -= sum * W;
1776 
1777  // CALCULATE THE JACOBIAN
1778  if (flag)
1779  {
1780  // Loop over the velocity shape functions again
1781  for (unsigned l2 = 0; l2 < n_node; l2++)
1782  {
1783  // Cannot loop over the velocity components
1784  // --- need to compute these separately instead
1785 
1786  unsigned i2 = 0;
1787 
1788  // If at a non-zero degree of freedom add in the entry
1789  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
1790  if (local_unknown >= 0)
1791  {
1792  // Convective r-term contribution
1793  double jac_conv = r * (psif(l2) * interpolated_dudx(0, 0) +
1794  u_r * dpsifdx(l2, 0));
1795 
1796  // Convective theta-term contribution
1797  jac_conv += u_theta * dpsifdx(l2, 1);
1798 
1799  // Add the time-derivative contribution and the convective
1800  // contribution to the sum
1801  double jac_sum =
1802  (scaled_re_st * node_pt(l2)->time_stepper_pt()->weight(1, 0) *
1803  psif(l2) * r2 +
1804  scaled_re * jac_conv * r) *
1805  testf(l);
1806 
1807  // Subtract the mesh-velocity terms
1808  if (!ALE_is_disabled)
1809  {
1810  jac_sum -= scaled_re_st *
1811  (mesh_velocity[0] * dpsifdx(l2, 0) * r +
1812  mesh_velocity[1] * dpsifdx(l2, 1)) *
1813  r * sin_theta * testf(l);
1814  }
1815 
1816 
1817  // Contribution from the r-derivative test function part
1818  // of stress tensor
1819  jac_sum += 2.0 * dpsifdx(l2, 0) * dtestfdx(l, 0) * r2;
1820 
1821  // Contribution from the theta-derivative test function part
1822  // of stress tensor
1823  jac_sum += dpsifdx(l2, 1) * dtestfdx(l, 1);
1824 
1825 
1826  // Contribution from the undifferentiated test function part
1827  // of stress tensor
1828  jac_sum += 4.0 * psif[l2] * testf(l);
1829 
1830  // Add the total contribution to the jacobian multiplied
1831  // by the jacobian weight
1832  jacobian(local_eqn, local_unknown) -= jac_sum * sin_theta * W;
1833 
1834  // Mass matrix term
1835  if (flag == 2)
1836  {
1837  mass_matrix(local_eqn, local_unknown) +=
1838  scaled_re_st * psif[l2] * testf[l] * r2 * sin_theta * W;
1839  }
1840  }
1841  // End of i2=0 section
1842 
1843 
1844  i2 = 1;
1845 
1846  // If at a non-zero degree of freedom add in the entry
1847  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
1848  if (local_unknown >= 0)
1849  {
1850  // No time derivative contribution
1851 
1852  // Convective contribution
1853  double jac_sum = scaled_re *
1854  (interpolated_dudx(0, 1) - 2.0 * u_theta) *
1855  psif(l2) * r * sin_theta * testf(l);
1856 
1857  // Mesh velocity terms
1858  if (!ALE_is_disabled)
1859  {
1860  jac_sum += scaled_re_st * mesh_velocity[1] * psif(l2) * r *
1861  sin_theta * testf(l);
1862  }
1863 
1864  // Contribution from the theta-derivative test function
1865  // part of stress tensor
1866  jac_sum +=
1867  (r * dpsifdx(l2, 0) - psif(l2)) * dtestfdx(l, 1) * sin_theta;
1868 
1869  // Contribution from the undifferentiated test function
1870  // part of stress tensor
1871  jac_sum += 2.0 *
1872  (dpsifdx(l2, 1) * sin_theta + psif(l2) * cos_theta) *
1873  testf(l);
1874 
1875  // Add the full contribution to the jacobian matrix
1876  jacobian(local_eqn, local_unknown) -= jac_sum * W;
1877 
1878  } // End of i2=1 section
1879 
1880 
1881  i2 = 2;
1882 
1883  // If at a non-zero degree of freedom add in the entry
1884  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
1885  if (local_unknown >= 0)
1886  {
1887  // Single convective-term contribution
1888  jacobian(local_eqn, local_unknown) += 2.0 * scaled_re * u_phi *
1889  psif(l2) * r * sin_theta *
1890  testf[l] * W;
1891 
1892  // Coriolis term
1893  jacobian(local_eqn, local_unknown) +=
1894  2.0 * scaled_re_inv_ro * sin_theta * psif(l2) * r2 *
1895  sin_theta * testf[l] * W;
1896  }
1897  // End of i2=2 section
1898  }
1899  // End of the l2 loop over the test functions
1900 
1901  // Now loop over pressure shape functions
1902  // This is the contribution from pressure gradient
1903  for (unsigned l2 = 0; l2 < n_pres; l2++)
1904  {
1905  /*If we are at a non-zero degree of freedom in the entry*/
1906  local_unknown = p_local_eqn(l2);
1907  if (local_unknown >= 0)
1908  {
1909  jacobian(local_eqn, local_unknown) +=
1910  psip(l2) * (r2 * dtestfdx(l, 0) + 2.0 * r * testf(l)) *
1911  sin_theta * W;
1912  }
1913  }
1914  } // End of Jacobian calculation
1915 
1916  } // End of if not boundary condition statement
1917  // End of r-component momentum equation
1918 
1919 
1920  // SECOND: theta-component momentum equation
1921 
1922  i = 1;
1923 
1924  // IF it's not a boundary condition
1925  local_eqn = nodal_local_eqn(l, u_nodal_index[i]);
1926  if (local_eqn >= 0)
1927  {
1928  // All convective terms
1929  double conv = (u_r * interpolated_dudx(1, 0) * r +
1930  u_theta * interpolated_dudx(1, 1) + u_r * u_theta) *
1931  sin_theta -
1932  u_phi * u_phi * cos_theta;
1933 
1934  // Add the time-derivative and convective terms to the
1935  // residuals sum
1936  double sum =
1937  (scaled_re_st * dudt[1] * r2 * sin_theta + scaled_re * r * conv) *
1938  testf(l);
1939 
1940  // Subtract the mesh velocity terms
1941  if (!ALE_is_disabled)
1942  {
1943  sum -= scaled_re_st *
1944  (mesh_velocity[0] * interpolated_dudx(1, 0) * r +
1945  mesh_velocity[1] * (interpolated_dudx(1, 1) + u_r)) *
1946  r * sin_theta * testf(l);
1947  }
1948 
1949  // Add the Coriolis term
1950  sum -= 2.0 * scaled_re_inv_ro * cos_theta * u_phi * r2 * sin_theta *
1951  testf[l];
1952 
1953  // r-derivative test function component of stress tensor
1954  sum +=
1955  (r * interpolated_dudx(1, 0) - u_theta + interpolated_dudx(0, 1)) *
1956  r * sin_theta * dtestfdx(l, 0);
1957 
1958  // theta-derivative test function component of stress tensor
1959  sum +=
1960  (-r * interpolated_p + 2.0 * interpolated_dudx(1, 1) + 2.0 * u_r) *
1961  dtestfdx(l, 1) * sin_theta;
1962 
1963  // Undifferentiated test function component of stress tensor
1964  sum -=
1965  ((r * interpolated_dudx(1, 0) - u_theta + interpolated_dudx(0, 1)) *
1966  sin_theta -
1967  (-r * interpolated_p + 2.0 * u_r + 2.0 * u_theta * cot_theta) *
1968  cos_theta) *
1969  testf(l);
1970 
1971  // Add the sum to the residuals,
1972  //(Negative for consistency)
1973  residuals[local_eqn] -= sum * W;
1974 
1975  // CALCULATE THE JACOBIAN
1976  if (flag)
1977  {
1978  // Loop over the velocity shape functions again
1979  for (unsigned l2 = 0; l2 < n_node; l2++)
1980  {
1981  // Cannot loop over the velocity components
1982  // --- need to compute these separately instead
1983 
1984  unsigned i2 = 0;
1985 
1986  // If at a non-zero degree of freedom add in the entry
1987  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
1988  if (local_unknown >= 0)
1989  {
1990  // No time derivative contribution
1991 
1992  // Convective terms
1993  double jac_sum = scaled_re *
1994  (r2 * interpolated_dudx(1, 0) + r * u_theta) *
1995  psif(l2) * sin_theta * testf(l);
1996 
1997  // Mesh velocity terms
1998  if (!ALE_is_disabled)
1999  {
2000  jac_sum -= scaled_re_st * mesh_velocity[1] * psif(l2) * r *
2001  sin_theta * testf(l);
2002  }
2003 
2004  // Contribution from the r-derivative
2005  // test function part of stress tensor
2006  jac_sum += dpsifdx(l2, 1) * dtestfdx(l, 0) * sin_theta * r;
2007 
2008  // Contribution from the theta-derivative test function
2009  // part of stress tensor
2010  jac_sum += 2.0 * psif(l2) * dtestfdx(l, 1) * sin_theta;
2011 
2012  // Contribution from the undifferentiated test function
2013  // part of stress tensor
2014  jac_sum -=
2015  (dpsifdx(l2, 1) * sin_theta - 2.0 * psif(l2) * cos_theta) *
2016  testf(l);
2017 
2018  // Add the sum to the jacobian
2019  jacobian(local_eqn, local_unknown) -= jac_sum * W;
2020  }
2021  // End of i2=0 section
2022 
2023 
2024  i2 = 1;
2025 
2026  // If at a non-zero degree of freedom add in the entry
2027  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
2028  if (local_unknown >= 0)
2029  {
2030  // Contribution from the convective terms
2031  double jac_conv = r * u_r * dpsifdx(l2, 0) +
2032  u_theta * dpsifdx(l2, 1) +
2033  (interpolated_dudx(1, 1) + u_r) * psif(l2);
2034 
2035  // Add the time-derivative term and the
2036  // convective terms
2037  double jac_sum =
2038  (scaled_re_st * node_pt(l2)->time_stepper_pt()->weight(1, 0) *
2039  psif(l2) * r2 +
2040  scaled_re * r * jac_conv) *
2041  testf(l) * sin_theta;
2042 
2043 
2044  // Mesh velocity terms
2045  if (!ALE_is_disabled)
2046  {
2047  jac_sum -= scaled_re_st *
2048  (mesh_velocity[0] * dpsifdx(l2, 0) * r +
2049  mesh_velocity[1] * dpsifdx(l2, 1)) *
2050  r * sin_theta * testf(l);
2051  }
2052 
2053  // Contribution from the r-derivative test function
2054  // part of stress tensor
2055  jac_sum += (r * dpsifdx(l2, 0) - psif(l2)) * r *
2056  dtestfdx(l, 0) * sin_theta;
2057 
2058  // Contribution from the theta-derivative test function
2059  // part of stress tensor
2060  jac_sum += 2.0 * dpsifdx(l2, 1) * dtestfdx(l, 1) * sin_theta;
2061 
2062  // Contribution from the undifferentiated test function
2063  // part of stress tensor
2064  jac_sum -= ((r * dpsifdx(l2, 0) - psif(l2)) * sin_theta -
2065  2.0 * psif(l2) * cot_theta * cos_theta) *
2066  testf(l);
2067 
2068  // Add the contribution to the jacobian
2069  jacobian(local_eqn, local_unknown) -= jac_sum * W;
2070 
2071  // Mass matrix term
2072  if (flag == 2)
2073  {
2074  mass_matrix(local_eqn, local_unknown) +=
2075  scaled_re_st * psif[l2] * testf[l] * r2 * sin_theta * W;
2076  }
2077  }
2078  // End of i2=1 section
2079 
2080 
2081  i2 = 2;
2082 
2083  // If at a non-zero degree of freedom add in the entry
2084  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
2085  if (local_unknown >= 0)
2086  {
2087  // Only a convective term contribution
2088  jacobian(local_eqn, local_unknown) += scaled_re * 2.0 * r *
2089  psif(l2) * u_phi *
2090  cos_theta * testf(l) * W;
2091 
2092  // Coriolis term
2093  jacobian(local_eqn, local_unknown) +=
2094  2.0 * scaled_re_inv_ro * cos_theta * psif(l2) * r2 *
2095  sin_theta * testf[l] * W;
2096  }
2097  // End of i2=2 section
2098  }
2099  // End of the l2 loop over the test functions
2100 
2101  /*Now loop over pressure shape functions*/
2102  /*This is the contribution from pressure gradient*/
2103  for (unsigned l2 = 0; l2 < n_pres; l2++)
2104  {
2105  /*If we are at a non-zero degree of freedom in the entry*/
2106  local_unknown = p_local_eqn(l2);
2107  if (local_unknown >= 0)
2108  {
2109  jacobian(local_eqn, local_unknown) +=
2110  psip(l2) * r *
2111  (dtestfdx(l, 1) * sin_theta + cos_theta * testf[l]) * W;
2112  }
2113  }
2114  } /*End of Jacobian calculation*/
2115 
2116  } // End of if not boundary condition statement
2117  // End of theta-component of momentum
2118 
2119 
2120  // THIRD: phi-component momentum equation
2121 
2122  i = 2;
2123 
2124  // IF it's not a boundary condition
2125  local_eqn = nodal_local_eqn(l, u_nodal_index[i]);
2126  if (local_eqn >= 0)
2127  {
2128  // Convective r-terms
2129  double conv = u_r * interpolated_dudx(2, 0) * r * sin_theta;
2130 
2131  // Convective theta-terms
2132  conv += u_theta * interpolated_dudx(2, 1) * sin_theta;
2133 
2134  // Remaining convective terms
2135  conv += u_phi * (u_r * sin_theta + u_theta * cos_theta);
2136 
2137  // Add the time-derivative term and the convective terms to the sum
2138  double sum =
2139  (scaled_re_st * r2 * dudt[2] * sin_theta + scaled_re * conv * r) *
2140  testf(l);
2141 
2142  // Mesh velocity terms
2143  if (!ALE_is_disabled)
2144  {
2145  sum -= scaled_re_st *
2146  (mesh_velocity[0] * interpolated_dudx(2, 0) * r +
2147  mesh_velocity[1] * interpolated_dudx(2, 1)) *
2148  r * sin_theta * testf(l);
2149  }
2150 
2151  // Add the Coriolis term
2152  sum += 2.0 * scaled_re_inv_ro *
2153  (cos_theta * u_theta + sin_theta * u_r) * r2 * sin_theta *
2154  testf[l];
2155 
2156  // r-derivative test function component of stress tensor
2157  sum += (r2 * interpolated_dudx(2, 0) - r * u_phi) * dtestfdx(l, 0) *
2158  sin_theta;
2159 
2160  // theta-derivative test function component of stress tensor
2161  sum += (interpolated_dudx(2, 1) * sin_theta - u_phi * cos_theta) *
2162  dtestfdx(l, 1);
2163 
2164  // Undifferentiated test function component of stress tensor
2165  sum -= ((r * interpolated_dudx(2, 0) - u_phi) * sin_theta +
2166  (interpolated_dudx(2, 1) - u_phi * cot_theta) * cos_theta) *
2167  testf(l);
2168 
2169  // Add the sum to the residuals
2170  residuals[local_eqn] -= sum * W;
2171 
2172 
2173  // CALCULATE THE JACOBIAN
2174  if (flag)
2175  {
2176  // Loop over the velocity shape functions again
2177  for (unsigned l2 = 0; l2 < n_node; l2++)
2178  {
2179  // Cannot loop over the velocity components
2180  // -- need to compute these separately instead
2181 
2182  unsigned i2 = 0;
2183 
2184  // If at a non-zero degree of freedom add in the entry
2185  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
2186  if (local_unknown >= 0)
2187  {
2188  // Contribution from convective terms
2189  jacobian(local_eqn, local_unknown) -=
2190  scaled_re * (r * interpolated_dudx(2, 0) + u_phi) * psif(l2) *
2191  testf(l) * r * sin_theta * W;
2192 
2193  // Coriolis term
2194  jacobian(local_eqn, local_unknown) -=
2195  2.0 * scaled_re_inv_ro * sin_theta * psif(l2) * r2 *
2196  sin_theta * testf[l] * W;
2197  }
2198  // End of i2=0 section
2199 
2200  i2 = 1;
2201 
2202  // If at a non-zero degree of freedom add in the entry
2203  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
2204  if (local_unknown >= 0)
2205  {
2206  // Contribution from convective terms
2207  jacobian(local_eqn, local_unknown) -=
2208  scaled_re *
2209  (interpolated_dudx(2, 1) * sin_theta + u_phi * cos_theta) *
2210  r * psif(l2) * testf(l) * W;
2211 
2212 
2213  // Coriolis term
2214  jacobian(local_eqn, local_unknown) -=
2215  2.0 * scaled_re_inv_ro * cos_theta * psif(l2) * r2 *
2216  sin_theta * testf[l] * W;
2217  }
2218  // End of i2=1 section
2219 
2220 
2221  i2 = 2;
2222 
2223  // If at a non-zero degree of freedom add in the entry
2224  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
2225  if (local_unknown >= 0)
2226  {
2227  // Convective terms
2228  double jac_conv = r * u_r * dpsifdx(l2, 0) * sin_theta;
2229 
2230  // Convective theta-term contribution
2231  jac_conv += u_theta * dpsifdx(l2, 1) * sin_theta;
2232 
2233  // Contribution from the remaining convective terms
2234  jac_conv += (u_r * sin_theta + u_theta * cos_theta) * psif(l2);
2235 
2236  // Add the convective and time derivatives
2237  double jac_sum =
2238  (scaled_re_st * node_pt(l2)->time_stepper_pt()->weight(1, 0) *
2239  psif(l2) * r2 * sin_theta +
2240  scaled_re * r * jac_conv) *
2241  testf(l);
2242 
2243 
2244  // Mesh velocity terms
2245  if (!ALE_is_disabled)
2246  {
2247  jac_sum -= scaled_re_st *
2248  (mesh_velocity[0] * dpsifdx(l2, 0) * r +
2249  mesh_velocity[1] * dpsifdx(l2, 1)) *
2250  r * sin_theta * testf(l);
2251  }
2252 
2253  // Contribution from the r-derivative test function
2254  // part of stress tensor
2255  jac_sum += (r * dpsifdx(l2, 0) - psif(l2)) * dtestfdx(l, 0) *
2256  r * sin_theta;
2257 
2258  // Contribution from the theta-derivative test function
2259  // part of stress tensor
2260  jac_sum += (dpsifdx(l2, 1) * sin_theta - psif(l2) * cos_theta) *
2261  dtestfdx(l, 1);
2262 
2263  // Contribution from the undifferentiated test function
2264  // part of stress tensor
2265  jac_sum -=
2266  ((r * dpsifdx(l2, 0) - psif(l2)) * sin_theta +
2267  (dpsifdx(l2, 1) - psif(l2) * cot_theta) * cos_theta) *
2268  testf(l);
2269 
2270  // Add to the jacobian
2271  jacobian(local_eqn, local_unknown) -= jac_sum * W;
2272 
2273  // Mass matrix term
2274  if (flag == 2)
2275  {
2276  mass_matrix(local_eqn, local_unknown) +=
2277  scaled_re_st * psif(l2) * testf[l] * r2 * sin_theta * W;
2278  }
2279  }
2280  // End of i2=2 section
2281  }
2282  // End of the l2 loop over the test functions
2283 
2284  // We assume phi-derivatives to be zero
2285  // No phi-contribution to the pressure gradient to include here
2286 
2287  } // End of Jacobian calculation
2288 
2289  } // End of if not boundary condition statement
2290  // End of phi-component of momentum
2291 
2292  } // End of loop over shape functions
2293 
2294 
2295  // CONTINUITY EQUATION
2296  //-------------------
2297 
2298  // Loop over the shape functions
2299  for (unsigned l = 0; l < n_pres; l++)
2300  {
2301  local_eqn = p_local_eqn(l);
2302  // If not a boundary conditions
2303  // Can't loop over the velocity components, so put these in separately
2304  if (local_eqn >= 0)
2305  {
2306  residuals[local_eqn] += ((2.0 * u_r + r * interpolated_dudx(0, 0) +
2307  interpolated_dudx(1, 1)) *
2308  sin_theta +
2309  u_theta * cos_theta) *
2310  r * testp(l) * W;
2311 
2312  // CALCULATE THE JACOBIAN
2313  if (flag)
2314  {
2315  // Loop over the velocity shape functions
2316  for (unsigned l2 = 0; l2 < n_node; l2++)
2317  {
2318  // Can't loop over velocity components, so put these in separately
2319 
2320  // Contribution from the r-component
2321  unsigned i2 = 0;
2322  /*If we're at a non-zero degree of freedom add it in*/
2323  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
2324 
2325  if (local_unknown >= 0)
2326  {
2327  jacobian(local_eqn, local_unknown) +=
2328  (2.0 * psif(l2) + r * dpsifdx(l2, 0)) * r * sin_theta *
2329  testp(l) * W;
2330  }
2331  // End of contribution from r-component
2332 
2333 
2334  // Contribution from the theta-component
2335  i2 = 1;
2336  /*If we're at a non-zero degree of freedom add it in*/
2337  local_unknown = nodal_local_eqn(l2, u_nodal_index[i2]);
2338 
2339  if (local_unknown >= 0)
2340  {
2341  jacobian(local_eqn, local_unknown) +=
2342  r * (dpsifdx(l2, 1) * sin_theta + psif(l2) * cos_theta) *
2343  testp(l) * W;
2344  }
2345  // End of contribution from theta-component
2346 
2347  } // End of loop over l2
2348  } // End of Jacobian calculation
2349 
2350  } // End of if not boundary condition
2351  } // End of loop over l
2352  }
2353  }
2354 
2355  /// ///////////////////////////////////////////////////////////////////
2356  /// ///////////////////////////////////////////////////////////////////
2357  /// ///////////////////////////////////////////////////////////////////
2358 
2359 
2360  //=========================================================================
2361  /// Add to the set \c paired_load_data pairs containing
2362  /// - the pointer to a Data object
2363  /// and
2364  /// - the index of the value in that Data object
2365  /// .
2366  /// for all values (pressures, velocities) that affect the
2367  /// load computed in the \c get_load(...) function.
2368  //=========================================================================
2370  std::set<std::pair<Data*, unsigned>>& paired_load_data)
2371  {
2372  // Find the index at which the velocity is stored
2373  unsigned u_index[3];
2374  for (unsigned i = 0; i < 3; i++)
2375  {
2376  u_index[i] = this->u_index_spherical_nst(i);
2377  }
2378 
2379  // Loop over the nodes
2380  unsigned n_node = this->nnode();
2381  for (unsigned n = 0; n < n_node; n++)
2382  {
2383  // Loop over the velocity components and add pointer to their data
2384  // and indices to the vectors
2385  for (unsigned i = 0; i < 3; i++)
2386  {
2387  paired_load_data.insert(std::make_pair(this->node_pt(n), u_index[i]));
2388  }
2389  }
2390 
2391  // Identify the pressure data
2392  this->identify_pressure_data(paired_load_data);
2393  }
2394 
2395  //=========================================================================
2396  /// Add to the set \c paired_pressure_data pairs containing
2397  /// - the pointer to a Data object
2398  /// and
2399  /// - the index of the value in that Data object
2400  /// .
2401  /// for pressures values that affect the
2402  /// load computed in the \c get_load(...) function.
2403  //=========================================================================
2405  std::set<std::pair<Data*, unsigned>>& paired_pressure_data)
2406  {
2407  // Loop over the internal data
2408  unsigned n_internal = this->ninternal_data();
2409  for (unsigned l = 0; l < n_internal; l++)
2410  {
2411  unsigned nval = this->internal_data_pt(l)->nvalue();
2412  // Add internal data
2413  for (unsigned j = 0; j < nval; j++)
2414  {
2415  paired_pressure_data.insert(
2416  std::make_pair(this->internal_data_pt(l), j));
2417  }
2418  }
2419  }
2420 
2421  //=============================================================================
2422  /// Create a list of pairs for all unknowns in this element,
2423  /// so that the first entry in each pair contains the global equation
2424  /// number of the unknown, while the second one contains the number
2425  /// of the "DOF type" that this unknown is associated with.
2426  /// (Function can obviously only be called if the equation numbering
2427  /// scheme has been set up.)
2428  //=============================================================================
2430  std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const
2431  {
2432  // number of nodes
2433  unsigned n_node = this->nnode();
2434 
2435  // number of pressure values
2436  unsigned n_press = this->npres_spherical_nst();
2437 
2438  // temporary pair (used to store dof lookup prior to being added to list)
2439  std::pair<unsigned, unsigned> dof_lookup;
2440 
2441  // pressure dof number
2442  unsigned pressure_dof_number = 3;
2443 
2444  // loop over the pressure values
2445  for (unsigned n = 0; n < n_press; n++)
2446  {
2447  // determine local eqn number
2448  int local_eqn_number = this->p_local_eqn(n);
2449 
2450  // ignore pinned values - far away degrees of freedom resulting
2451  // from hanging nodes can be ignored since these are be dealt
2452  // with by the element containing their master nodes
2453  if (local_eqn_number >= 0)
2454  {
2455  // store dof lookup in temporary pair: First entry in pair
2456  // is global equation number; second entry is dof type
2457  dof_lookup.first = this->eqn_number(local_eqn_number);
2458  dof_lookup.second = pressure_dof_number;
2459 
2460  // add to list
2461  dof_lookup_list.push_front(dof_lookup);
2462  }
2463  }
2464 
2465  // loop over the nodes
2466  for (unsigned n = 0; n < n_node; n++)
2467  {
2468  // find the number of values at this node
2469  unsigned nv = this->node_pt(n)->nvalue();
2470 
2471  // loop over these values
2472  for (unsigned v = 0; v < nv; v++)
2473  {
2474  // determine local eqn number
2475  int local_eqn_number = this->nodal_local_eqn(n, v);
2476 
2477  // ignore pinned values
2478  if (local_eqn_number >= 0)
2479  {
2480  // store dof lookup in temporary pair: First entry in pair
2481  // is global equation number; second entry is dof type
2482  dof_lookup.first = this->eqn_number(local_eqn_number);
2483  dof_lookup.second = v;
2484 
2485  // add to list
2486  dof_lookup_list.push_front(dof_lookup);
2487  }
2488  }
2489  }
2490  }
2491 
2492 
2493  /// ////////////////////////////////////////////////////////////////////////
2494  /// ////////////////////////////////////////////////////////////////////////
2495  /// ////////////////////////////////////////////////////////////////////////
2496 
2497  // 2D Taylor--Hood
2498  // Set the data for the number of Variables at each node
2499  const unsigned QSphericalTaylorHoodElement::Initial_Nvalue[9] = {
2500  4, 3, 4, 3, 3, 3, 4, 3, 4};
2501 
2502  // Set the data for the pressure conversion array
2503  const unsigned QSphericalTaylorHoodElement::Pconv[4] = {0, 2, 6, 8};
2504 
2505  //=========================================================================
2506  /// Add to the set \c paired_load_data pairs containing
2507  /// - the pointer to a Data object
2508  /// and
2509  /// - the index of the value in that Data object
2510  /// .
2511  /// for all values (pressures, velocities) that affect the
2512  /// load computed in the \c get_load(...) function.
2513  //=========================================================================
2515  std::set<std::pair<Data*, unsigned>>& paired_load_data)
2516  {
2517  // Find the index at which the velocity is stored
2518  unsigned u_index[3];
2519  for (unsigned i = 0; i < 3; i++)
2520  {
2521  u_index[i] = this->u_index_spherical_nst(i);
2522  }
2523 
2524  // Loop over the nodes
2525  unsigned n_node = this->nnode();
2526  for (unsigned n = 0; n < n_node; n++)
2527  {
2528  // Loop over the velocity components and add pointer to their data
2529  // and indices to the vectors
2530  for (unsigned i = 0; i < 3; i++)
2531  {
2532  paired_load_data.insert(std::make_pair(this->node_pt(n), u_index[i]));
2533  }
2534  }
2535 
2536  // Identify the pressure data
2537  this->identify_pressure_data(paired_load_data);
2538  }
2539 
2540  //=========================================================================
2541  /// Add to the set \c paired_load_data pairs containing
2542  /// - the pointer to a Data object
2543  /// and
2544  /// - the index of the value in that Data object
2545  /// .
2546  /// for pressure values that affect the
2547  /// load computed in the \c get_load(...) function.
2548  //=========================================================================
2550  std::set<std::pair<Data*, unsigned>>& paired_pressure_data)
2551  {
2552  // Find the index at which the pressure is stored
2553  unsigned p_index =
2554  static_cast<unsigned>(this->p_nodal_index_spherical_nst());
2555 
2556  // Loop over the pressure data
2557  unsigned n_pres = npres_spherical_nst();
2558  for (unsigned l = 0; l < n_pres; l++)
2559  {
2560  // The DIMth entry in each nodal data is the pressure, which
2561  // affects the traction
2562  paired_pressure_data.insert(
2563  std::make_pair(this->node_pt(Pconv[l]), p_index));
2564  }
2565  }
2566 
2567 
2568  //============================================================================
2569  /// Create a list of pairs for all unknowns in this element,
2570  /// so the first entry in each pair contains the global equation
2571  /// number of the unknown, while the second one contains the number
2572  /// of the "DOF type" that this unknown is associated with.
2573  /// (Function can obviously only be called if the equation numbering
2574  /// scheme has been set up.)
2575  //============================================================================
2577  std::list<std::pair<unsigned long, unsigned>>& dof_lookup_list) const
2578  {
2579  // number of nodes
2580  unsigned n_node = this->nnode();
2581 
2582  // local eqn no for pressure unknown
2583  // unsigned p_index = this->p_nodal_index_spherical_nst();
2584 
2585  // temporary pair (used to store dof lookup prior to being added to list)
2586  std::pair<unsigned, unsigned> dof_lookup;
2587 
2588  // loop over the nodes
2589  for (unsigned n = 0; n < n_node; n++)
2590  {
2591  // find the number of values at this node
2592  unsigned nv = this->node_pt(n)->nvalue();
2593 
2594  // loop over these values
2595  for (unsigned v = 0; v < nv; v++)
2596  {
2597  // determine local eqn number
2598  int local_eqn_number = this->nodal_local_eqn(n, v);
2599 
2600  // ignore pinned values - far away degrees of freedom resulting
2601  // from hanging nodes can be ignored since these are be dealt
2602  // with by the element containing their master nodes
2603  if (local_eqn_number >= 0)
2604  {
2605  // store dof lookup in temporary pair: Global equation number
2606  // is the first entry in pair
2607  dof_lookup.first = this->eqn_number(local_eqn_number);
2608 
2609  // set dof numbers: Dof number is the second entry in pair
2610  dof_lookup.second = v;
2611 
2612  // add to list
2613  dof_lookup_list.push_front(dof_lookup);
2614  }
2615  }
2616  }
2617  }
2618 
2619 
2620 } // namespace oomph
static char t char * s
Definition: cfortran.h:568
cstr elem_len * i
Definition: cfortran.h:603
char t
Definition: cfortran.h:568
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
unsigned nvalue() const
Return number of values stored in data object (incl pinned ones).
Definition: nodes.h:483
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:4103
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2175
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:2593
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:3161
virtual double dshape_eulerian_at_knot(const unsigned &ipt, Shape &psi, DShape &dpsidx) const
Return the geometric shape functions and also first derivatives w.r.t. global coordinates at the ipt-...
Definition: elements.cc:3325
virtual void shape(const Vector< double > &s, Shape &psi) const =0
Calculate the geometric shape functions at local coordinate s. This function must be overloaded for e...
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:3962
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:1432
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2210
void(* SteadyExactSolutionFctPt)(const Vector< double > &, Vector< double > &)
Function pointer for function that computes vector-valued steady "exact solution" as .
Definition: elements.h:1759
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1963
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:3148
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:3186
double nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. If the node is hanging, the appropriate interpolation is ...
Definition: elements.h:2317
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:3298
virtual double J_eulerian_at_knot(const unsigned &ipt) const
Return the Jacobian of the mapping from local to global coordinates at the ipt-th integration point.
Definition: elements.cc:4168
double raw_dnodal_position_dt(const unsigned &n, const unsigned &i) const
Return the i-th component of nodal velocity: dx/dt at local node n. Do not use the hanging node repre...
Definition: elements.h:2256
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:2576
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:1686
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:3174
void(* UnsteadyExactSolutionFctPt)(const double &, const Vector< double > &, Vector< double > &)
Function pointer for function that computes Vector-valued time-dependent function as .
Definition: elements.h:1765
double dnodal_position_dt(const unsigned &n, const unsigned &i) const
Return the i-th component of nodal velocity: dx/dt at local node n.
Definition: elements.h:2333
virtual void shape_at_knot(const unsigned &ipt, Shape &psi) const
Return the geometric shape function at the ipt-th integration point.
Definition: elements.cc:3220
unsigned ndof() const
Return the number of equations/dofs in the element.
Definition: elements.h:835
unsigned long eqn_number(const unsigned &ieqn_local) const
Return the global equation number corresponding to the ieqn_local-th local equation number.
Definition: elements.h:704
int local_eqn_number(const unsigned long &ieqn_global) const
Return the local equation number corresponding to the ieqn_global-th global equation number....
Definition: elements.h:726
Data *& internal_data_pt(const unsigned &i)
Return a pointer to i-th internal data object.
Definition: elements.h:622
unsigned ninternal_data() const
Return the number of internal data objects.
Definition: elements.h:823
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.
An OomphLibError object which should be thrown when an run-time error is encountered....
An OomphLibWarning object which should be created as a temporary object to issue a warning....
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Add to the set paired_load_data pairs containing.
int p_local_eqn(const unsigned &n) const
Return the local equation numbers for the pressure values.
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Add to the set paired_pressure_data pairs containing.
unsigned npres_spherical_nst() const
Return number of pressure values.
int p_nodal_index_spherical_nst() const
Set the value at which the pressure is stored in the nodes In this case the third index because there...
void get_dof_numbers_for_unknowns(std::list< std::pair< unsigned long, unsigned >> &dof_lookup_list) const
Create a list of pairs for all unknowns in this element, so that the first entry in each pair contain...
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
void identify_load_data(std::set< std::pair< Data *, unsigned >> &paired_load_data)
Add to the set paired_load_data pairs containing.
static const unsigned Pconv[]
Static array of ints to hold conversion from pressure node numbers to actual node numbers.
void identify_pressure_data(std::set< std::pair< Data *, unsigned >> &paired_pressure_data)
Add to the set paired_pressure_data pairs containing.
unsigned npres_spherical_nst() const
Return number of pressure values.
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
double d_kin_energy_dt() const
Get integral of time derivative of kinetic energy over element.
void compute_error(std::ostream &outfile, FiniteElement::UnsteadyExactSolutionFctPt exact_soln_pt, const double &time, double &error, double &norm)
Validate against exact solution at given time Solution is provided via function pointer....
void get_vorticity(const Vector< double > &s, Vector< double > &vorticity) const
Compute the vorticity vector at local coordinate s.
double interpolated_dudx_spherical_nst(const Vector< double > &s, const unsigned &i, const unsigned &j) const
Return matrix entry dudx(i,j) of the FE interpolated velocity derivative at local coordinate s.
virtual double dshape_and_dtest_eulerian_at_knot_spherical_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and derivatives w.r.t. global coords at ipt-th integration point Return J...
double pressure_integral() const
Integral of pressure over element.
virtual void pshape_spherical_nst(const Vector< double > &s, Shape &psi) const =0
Compute the pressure shape functions at local coordinate s.
void compute_error_e(std::ostream &outfile, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt, FiniteElement::SteadyExactSolutionFctPt exact_soln_dr_pt, FiniteElement::SteadyExactSolutionFctPt exact_soln_dtheta_pt, double &u_error, double &u_norm, double &p_error, double &p_norm)
Validate against exact solution. Solution is provided direct from exact_soln function....
void output(std::ostream &outfile)
Output function: x,y,[z],u,v,[w],p in tecplot format. Default number of plot points.
void full_output(std::ostream &outfile)
Full output function: x,y,[z],u,v,[w],p,du/dt,dv/dt,[dw/dt],dissipation in tecplot format....
virtual unsigned u_index_spherical_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
const double & re() const
Reynolds number.
static double Default_Physical_Ratio_Value
Static default value for the physical ratios (all are initialised to one)
void get_pressure_and_velocity_mass_matrix_diagonal(Vector< double > &press_mass_diag, Vector< double > &veloc_mass_diag, const unsigned &which_one=0)
Compute the diagonal of the velocity/pressure mass matrices. If which one=0, both are computed,...
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Strain-rate tensor: 1/2 (du_i/dx_j + du_j/dx_i)
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
const double & re_invro() const
Global Reynolds number multiplied by inverse Rossby number.
virtual double p_spherical_nst(const unsigned &n_p) const =0
Pressure at local pressure "node" n_p Uses suitably interpolated value for hanging nodes.
const double & density_ratio() const
Density ratio for element: Element's density relative to the viscosity used in the definition of the ...
virtual void get_body_force_spherical_nst(const double &time, const unsigned &ipt, const Vector< double > &s, const Vector< double > &x, Vector< double > &result)
Calculate the body force at a given time and local and/or Eulerian position. This function is virtual...
static int Pressure_not_stored_at_node
Static "magic" number that indicates that the pressure is not stored at a node.
double du_dt_spherical_nst(const unsigned &n, const unsigned &i) const
i-th component of du/dt at local node n. Uses suitably interpolated value for hanging nodes.
double cot(const double &th) const
Include a cot function to simplify the NS momentum and jacobian expressions.
double interpolated_p_spherical_nst(const Vector< double > &s) const
Return FE interpolated pressure at local coordinate s.
void output_fct(std::ostream &outfile, const unsigned &nplot, FiniteElement::SteadyExactSolutionFctPt exact_soln_pt)
Output exact solution specified via function pointer at a given number of plot points....
double kin_energy() const
Get integral of kinetic energy over element.
void output_vorticity(std::ostream &outfile, const unsigned &nplot)
Output function: x,y,[z], [omega_x,omega_y,[and/or omega_z]] in tecplot format. nplot points in each ...
double dissipation() const
Return integral of dissipation over element.
static double Default_Physical_Constant_Value
Static default value for the physical constants (all initialised to zero)
void get_traction(const Vector< double > &s, const Vector< double > &N, Vector< double > &traction)
Compute traction (on the viscous scale) exerted onto the fluid at local coordinate s....
const Vector< double > & g() const
Vector of gravitational components.
virtual void fill_in_generic_residual_contribution_spherical_nst(Vector< double > &residuals, DenseMatrix< double > &jacobian, DenseMatrix< double > &mass_matrix, unsigned flag)
Compute the residuals for the Navier–Stokes equations; flag=1(or 0): do (or don't) compute the Jacobi...
static Vector< double > Default_Gravity_vector
Static default value for the gravity vector.
virtual int p_local_eqn(const unsigned &n) const =0
Access function for the local equation number information for the pressure. p_local_eqn[n] = local eq...
virtual unsigned npres_spherical_nst() const =0
Function to return number of pressure degrees of freedom.
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when time-derivatives are computed....
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Output function: x,y,[z],u,v,[w] in tecplot format. nplot points in each coordinate direction at time...
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
void interpolated_u_spherical_nst(const Vector< double > &s, Vector< double > &veloc) const
Compute vector of FE interpolated velocity u at local coordinate s.
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
Time *const & time_pt() const
Access function for the pointer to time (const version)
Definition: timesteppers.h:572
double & time()
Return the current value of the continuous time.
Definition: timesteppers.h:123
const double Pi
50 digits from maple
//////////////////////////////////////////////////////////////////// ////////////////////////////////...