linearised_axisym_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-2024 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Non-inline functions for linearised axisymmetric Navier-Stokes elements
27 
28 // oomph-lib includes
30 
31 namespace oomph
32 {
33  //=======================================================================
34  /// Linearised axisymmetric Navier--Stokes equations static data
35  //=======================================================================
36 
37  // Use the stress-divergence form by default (Gamma=1)
39 
40  // "Magic" number to indicate pressure is not stored at node
42  -100;
43 
44  // Physical constants default to zero
45  double LinearisedAxisymmetricNavierStokesEquations::
46  Default_Physical_Constant_Value = 0.0;
47 
48  // Azimuthal mode number defaults to zero
49  int LinearisedAxisymmetricNavierStokesEquations::
50  Default_Azimuthal_Mode_Number_Value = 0;
51 
52  // Density/viscosity ratios default to one
53  double
55  1.0;
56 
57 
58  //=======================================================================
59  /// Output function in tecplot format: Velocities only
60  /// r, z, U^C, U^S, V^C, V^S, W^C, W^S
61  /// at specified previous timestep (t=0 present; t>0 previous timestep).
62  /// Specified number of plot points in each coordinate direction.
63  //=======================================================================
65  std::ostream& outfile, const unsigned& nplot, const unsigned& t)
66  {
67  // Determine number of nodes in element
68  const unsigned n_node = nnode();
69 
70  // Provide storage for local shape functions
71  Shape psi(n_node);
72 
73  // Provide storage for vectors of local coordinates and
74  // global coordinates and velocities
75  Vector<double> s(2);
77  Vector<double> interpolated_u(6);
78 
79  // Tecplot header info
80  outfile << tecplot_zone_string(nplot);
81 
82  // Determine number of plot points
83  const unsigned n_plot_points = nplot_points(nplot);
84 
85  // Loop over plot points
86  for (unsigned iplot = 0; iplot < n_plot_points; iplot++)
87  {
88  // Get local coordinates of plot point
89  get_s_plot(iplot, nplot, s);
90 
91  // Get shape functions
92  shape(s, psi);
93 
94  // Loop over coordinate directions
95  for (unsigned i = 0; i < 2; i++)
96  {
97  // Initialise global coordinate
98  interpolated_x[i] = 0.0;
99 
100  // Loop over the local nodes and sum
101  for (unsigned l = 0; l < n_node; l++)
102  {
103  interpolated_x[i] += nodal_position(t, l, i) * psi[l];
104  }
105  }
106 
107  // Loop over the velocity components
108  for (unsigned i = 0; i < 6; i++)
109  {
110  // Get the index at which the velocity is stored
111  const unsigned u_nodal_index = u_index_linearised_axi_nst(i);
112 
113  // Initialise velocity
114  interpolated_u[i] = 0.0;
115 
116  // Loop over the local nodes and sum
117  for (unsigned l = 0; l < n_node; l++)
118  {
119  interpolated_u[i] += nodal_value(t, l, u_nodal_index) * psi[l];
120  }
121  }
122 
123  // Output global coordinates to file
124  for (unsigned i = 0; i < 2; i++)
125  {
126  outfile << interpolated_x[i] << " ";
127  }
128 
129  // Output velocities to file
130  for (unsigned i = 0; i < 6; i++)
131  {
132  outfile << interpolated_u[i] << " ";
133  }
134 
135  outfile << std::endl;
136  }
137 
138  // Write tecplot footer (e.g. FE connectivity lists)
139  write_tecplot_zone_footer(outfile, nplot);
140 
141  } // End of output_veloc
142 
143 
144  //=======================================================================
145  /// Output function in tecplot format:
146  /// r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S
147  /// in tecplot format. Specified number of plot points in each
148  /// coordinate direction.
149  //=======================================================================
151  std::ostream& outfile, const unsigned& nplot)
152  {
153  // Provide storage for vector of local coordinates
154  Vector<double> s(2);
155 
156  // Tecplot header info
157  outfile << tecplot_zone_string(nplot);
158 
159  // Determine number of plot points
160  const unsigned n_plot_points = nplot_points(nplot);
161 
162  // Loop over plot points
163  for (unsigned iplot = 0; iplot < n_plot_points; iplot++)
164  {
165  // Get local coordinates of plot point
166  get_s_plot(iplot, nplot, s);
167 
168  // Output global coordinates to file
169  for (unsigned i = 0; i < 2; i++)
170  {
171  outfile << interpolated_x(s, i) << " ";
172  }
173 
174  // Output velocities to file
175  for (unsigned i = 0; i < 6; i++)
176  {
177  outfile << interpolated_u_linearised_axi_nst(s, i) << " ";
178  }
179 
180  // Output pressure to file
181  for (unsigned i = 0; i < 2; i++)
182  {
183  outfile << interpolated_p_linearised_axi_nst(s, i) << " ";
184  }
185 
186  outfile << std::endl;
187  }
188  outfile << std::endl;
189 
190  // Write tecplot footer (e.g. FE connectivity lists)
191  write_tecplot_zone_footer(outfile, nplot);
192 
193  } // End of output
194 
195 
196  //=======================================================================
197  /// Output function in tecplot format:
198  /// r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S
199  /// Specified number of plot points in each coordinate direction.
200  //=======================================================================
202  FILE* file_pt, const unsigned& nplot)
203  {
204  // Provide storage for vector of local coordinates
205  Vector<double> s(2);
206 
207  // Tecplot header info
208  fprintf(file_pt, "%s ", tecplot_zone_string(nplot).c_str());
209 
210  // Determine number of plot points
211  const unsigned n_plot_points = nplot_points(nplot);
212 
213  // Loop over plot points
214  for (unsigned iplot = 0; iplot < n_plot_points; iplot++)
215  {
216  // Get local coordinates of plot point
217  get_s_plot(iplot, nplot, s);
218 
219  // Output global coordinates to file
220  for (unsigned i = 0; i < 2; i++)
221  {
222  fprintf(file_pt, "%g ", interpolated_x(s, i));
223  }
224 
225  // Output velocities to file
226  for (unsigned i = 0; i < 6; i++)
227  {
228  fprintf(file_pt, "%g ", interpolated_u_linearised_axi_nst(s, i));
229  }
230 
231  // Output pressure to file
232  for (unsigned i = 0; i < 2; i++)
233  {
234  fprintf(file_pt, "%g ", interpolated_p_linearised_axi_nst(s, i));
235  }
236 
237  fprintf(file_pt, "\n");
238  }
239 
240  fprintf(file_pt, "\n");
241 
242  // Write tecplot footer (e.g. FE connectivity lists)
243  write_tecplot_zone_footer(file_pt, nplot);
244 
245  } // End of output
246 
247 
248  //=======================================================================
249  /// Get strain-rate tensor: \f$ e_{ij} \f$ where
250  /// \f$ i,j = r,z,\theta \f$ (in that order). We evaluate this tensor
251  /// at a value of theta such that the product of theta and the azimuthal
252  /// mode number (k) gives \f$ \pi/4 \f$. Therefore
253  /// \f$ \cos(k \theta) = \sin(k \theta) = 1/\sqrt{2} \f$.
254  //=======================================================================
256  const Vector<double>& s, DenseMatrix<double>& strainrate) const
257  {
258 #ifdef PARANOID
259  if ((strainrate.ncol() != 3) || (strainrate.nrow() != 3))
260  {
261  std::ostringstream error_message;
262  error_message << "The strain rate has incorrect dimensions "
263  << strainrate.ncol() << " x " << strainrate.nrow()
264  << " Not 3" << std::endl;
265 
266  throw OomphLibError(
267  error_message.str(),
268  "LinearisedAxisymmetricNavierStokeEquations::strain_rate()",
269  OOMPH_EXCEPTION_LOCATION);
270  }
271 #endif
272 
273  // Determine number of nodes in element
274  const unsigned n_node = nnode();
275 
276  // Set up memory for the shape and test functions
277  Shape psi(n_node);
278  DShape dpsidx(n_node, 2);
279 
280  // Call the derivatives of the shape functions
281  dshape_eulerian(s, psi, dpsidx);
282 
283  // Radius
284  double interpolated_r = 0.0;
285 
286  // Velocity components and their derivatives
287  double UC = 0.0, US = 0.0;
288  double dUCdr = 0.0, dUSdr = 0.0;
289  double dUCdz = 0.0, dUSdz = 0.0;
290  double WC = 0.0, WS = 0.0;
291  double dWCdr = 0.0, dWSdr = 0.0;
292  double dWCdz = 0.0, dWSdz = 0.0;
293  double VC = 0.0, VS = 0.0;
294  double dVCdr = 0.0, dVSdr = 0.0;
295  double dVCdz = 0.0, dVSdz = 0.0;
296 
297  // Get the local storage for the indices
298  unsigned u_nodal_index[6];
299  for (unsigned i = 0; i < 6; ++i)
300  {
301  u_nodal_index[i] = u_index_linearised_axi_nst(i);
302  }
303 
304  // Loop over nodes to assemble velocities and their derivatives
305  // w.r.t. r and z (x_0 and x_1)
306  for (unsigned l = 0; l < n_node; l++)
307  {
308  interpolated_r += nodal_position(l, 0) * psi[l];
309 
310  UC += nodal_value(l, u_nodal_index[0]) * psi[l];
311  US += nodal_value(l, u_nodal_index[1]) * psi[l];
312  WC += nodal_value(l, u_nodal_index[2]) * psi[l];
313  WS += nodal_value(l, u_nodal_index[3]) * psi[l];
314  VC += nodal_value(l, u_nodal_index[4]) * psi[l];
315  VS += nodal_value(l, u_nodal_index[4]) * psi[l];
316 
317  dUCdr += nodal_value(l, u_nodal_index[0]) * dpsidx(l, 0);
318  dUSdr += nodal_value(l, u_nodal_index[1]) * dpsidx(l, 0);
319  dWCdr += nodal_value(l, u_nodal_index[2]) * dpsidx(l, 0);
320  dWSdr += nodal_value(l, u_nodal_index[3]) * dpsidx(l, 0);
321  dVCdr += nodal_value(l, u_nodal_index[4]) * dpsidx(l, 0);
322  dVSdr += nodal_value(l, u_nodal_index[5]) * dpsidx(l, 0);
323 
324  dUCdz += nodal_value(l, u_nodal_index[0]) * dpsidx(l, 1);
325  dUSdz += nodal_value(l, u_nodal_index[1]) * dpsidx(l, 1);
326  dWCdz += nodal_value(l, u_nodal_index[2]) * dpsidx(l, 1);
327  dWSdz += nodal_value(l, u_nodal_index[3]) * dpsidx(l, 1);
328  dVCdz += nodal_value(l, u_nodal_index[4]) * dpsidx(l, 1);
329  dVSdz += nodal_value(l, u_nodal_index[5]) * dpsidx(l, 1);
330  }
331 
332  // Cache azimuthal mode number
333  const int k = this->azimuthal_mode_number();
334 
335  // We wish to evaluate the strain-rate tensor at a value of theta
336  // such that k*theta = pi/4 radians. That way we pick up equal
337  // contributions from the real and imaginary parts of the velocities.
338  // sin(pi/4) = cos(pi/4) = 1/sqrt(2)
339  const double cosktheta = 1.0 / sqrt(2);
340  const double sinktheta = cosktheta;
341 
342  // Assemble velocities and their derivatives w.r.t. r, z and theta
343  // from real and imaginary parts
344  const double ur = UC * cosktheta + US * sinktheta;
345  const double utheta = VC * cosktheta + VS * sinktheta;
346 
347  const double durdr = dUCdr * cosktheta + dUSdr * sinktheta;
348  const double durdz = dUCdz * cosktheta + dUSdz * sinktheta;
349  const double durdtheta = k * US * cosktheta - k * UC * sinktheta;
350 
351  const double duzdr = dWCdr * cosktheta + dWSdr * sinktheta;
352  const double duzdz = dWCdz * cosktheta + dWSdz * sinktheta;
353  const double duzdtheta = k * WS * cosktheta - k * WC * sinktheta;
354 
355  const double duthetadr = dVCdr * cosktheta + dVSdr * sinktheta;
356  const double duthetadz = dVCdz * cosktheta + dVSdz * sinktheta;
357  const double duthetadtheta = k * VS * cosktheta - k * VC * sinktheta;
358 
359  // Assign strain rates without negative powers of the radius
360  // and zero those with:
361  strainrate(0, 0) = durdr;
362  strainrate(0, 1) = 0.5 * (durdz + duzdr);
363  strainrate(1, 0) = strainrate(0, 1);
364  strainrate(0, 2) = 0.5 * duthetadr;
365  strainrate(2, 0) = strainrate(0, 2);
366  strainrate(1, 1) = duzdz;
367  strainrate(1, 2) = 0.5 * duthetadz;
368  strainrate(2, 1) = strainrate(1, 2);
369  strainrate(2, 2) = 0.0;
370 
371  // Overwrite the strain rates with negative powers of the radius
372  // unless we're at the origin
373  if (std::abs(interpolated_r) > 1.0e-16)
374  {
375  double inverse_radius = 1.0 / interpolated_r;
376  strainrate(0, 2) =
377  0.5 * (duthetadr + inverse_radius * (durdtheta - utheta));
378  strainrate(2, 0) = strainrate(0, 2);
379  strainrate(2, 2) = inverse_radius * (ur + duthetadtheta);
380  strainrate(1, 2) = 0.5 * (duthetadz + inverse_radius * duzdtheta);
381  strainrate(2, 1) = strainrate(1, 2);
382  }
383 
384  } // End of strain_rate
385 
386 
387  //=======================================================================
388  /// Compute the residuals for the linearised axisymmetric Navier--Stokes
389  /// equations; flag=1(or 0): do (or don't) compute the Jacobian as well.
390  //=======================================================================
393  Vector<double>& residuals,
394  DenseMatrix<double>& jacobian,
395  DenseMatrix<double>& mass_matrix,
396  unsigned flag)
397  {
398  // Get the time from the first node in the element
399  const double time = this->node_pt(0)->time_stepper_pt()->time();
400 
401  // Determine number of nodes in the element
402  const unsigned n_node = nnode();
403 
404  // Determine how many pressure values there are associated with
405  // a single pressure component
406  const unsigned n_pres = npres_linearised_axi_nst();
407 
408  // Get the nodal indices at which the velocity is stored
409  unsigned u_nodal_index[6];
410  for (unsigned i = 0; i < 6; ++i)
411  {
412  u_nodal_index[i] = u_index_linearised_axi_nst(i);
413  }
414 
415  // Set up memory for the fluid shape and test functions
416  // Note that there are two spatial dimensions, r and z, in this problem
417  Shape psif(n_node), testf(n_node);
418  DShape dpsifdx(n_node, 2), dtestfdx(n_node, 2);
419 
420  // Set up memory for the pressure shape and test functions
421  Shape psip(n_pres), testp(n_pres);
422 
423  // Determine number of integration points
424  const unsigned n_intpt = integral_pt()->nweight();
425 
426  // Set up memory for the vector to hold local coordinates (two dimensions)
427  Vector<double> s(2);
428 
429  // Get physical variables from the element
430  // (Reynolds number must be multiplied by the density ratio)
431  const double scaled_re = re() * density_ratio();
432  const double scaled_re_st = re_st() * density_ratio();
433  const double visc_ratio = viscosity_ratio();
434  const int k = azimuthal_mode_number();
435 
436  // Integers used to store the local equation and unknown numbers
437  int local_eqn = 0, local_unknown = 0;
438 
439  // Loop over the integration points
440  for (unsigned ipt = 0; ipt < n_intpt; ipt++)
441  {
442  // Assign values of the local coordinates s
443  for (unsigned i = 0; i < 2; i++)
444  {
445  s[i] = integral_pt()->knot(ipt, i);
446  }
447 
448  // Get the integral weight
449  const double w = integral_pt()->weight(ipt);
450 
451  // Calculate the fluid shape and test functions, and their derivatives
452  // w.r.t. the global coordinates
454  ipt, psif, dpsifdx, testf, dtestfdx);
455 
456  // Calculate the pressure shape and test functions
457  pshape_linearised_axi_nst(s, psip, testp);
458 
459  // Premultiply the weights and the Jacobian of the mapping between
460  // local and global coordinates
461  const double W = w * J;
462 
463  // Allocate storage for the position and the derivative of the
464  // mesh positions w.r.t. time
466  Vector<double> mesh_velocity(2, 0.0);
467 
468  // Allocate storage for the velocity components (six of these)
469  // and their derivatives w.r.t. time
470  Vector<double> interpolated_u(6, 0.0);
471  Vector<double> dudt(6, 0.0);
472 
473  // Allocate storage for the pressure components (two of these)
474  Vector<double> interpolated_p(2, 0.0);
475 
476  // Allocate storage for the derivatives of the velocity components
477  // w.r.t. global coordinates (r and z)
478  DenseMatrix<double> interpolated_dudx(6, 2, 0.0);
479 
480  // Calculate pressure at the integration point
481  // -------------------------------------------
482 
483  // Loop over pressure degrees of freedom (associated with a single
484  // pressure component) in the element
485  for (unsigned l = 0; l < n_pres; l++)
486  {
487  // Cache the shape function
488  const double psip_ = psip(l);
489 
490  // Loop over the two pressure components
491  for (unsigned i = 0; i < 2; i++)
492  {
493  // Get the value
494  const double p_value = this->p_linearised_axi_nst(l, i);
495 
496  // Add contribution
497  interpolated_p[i] += p_value * psip_;
498  }
499  } // End of loop over the pressure degrees of freedom in the element
500 
501  // Calculate velocities and their derivatives at the integration point
502  // -------------------------------------------------------------------
503 
504  // Loop over the element's nodes
505  for (unsigned l = 0; l < n_node; l++)
506  {
507  // Cache the shape function
508  const double psif_ = psif(l);
509 
510  // Loop over the two coordinate directions
511  for (unsigned i = 0; i < 2; i++)
512  {
513  interpolated_x[i] += this->raw_nodal_position(l, i) * psif_;
514  }
515 
516  // Loop over the six velocity components
517  for (unsigned i = 0; i < 6; i++)
518  {
519  // Get the value
520  const double u_value = this->raw_nodal_value(l, u_nodal_index[i]);
521 
522  // Add contribution
523  interpolated_u[i] += u_value * psif_;
524 
525  // Add contribution to dudt
526  dudt[i] += du_dt_linearised_axi_nst(l, i) * psif_;
527 
528  // Loop over two coordinate directions (for derivatives)
529  for (unsigned j = 0; j < 2; j++)
530  {
531  interpolated_dudx(i, j) += u_value * dpsifdx(l, j);
532  }
533  }
534  } // End of loop over the element's nodes
535 
536  // Get the mesh velocity if ALE is enabled
537  if (!ALE_is_disabled)
538  {
539  // Loop over the element's nodes
540  for (unsigned l = 0; l < n_node; l++)
541  {
542  // Loop over the two coordinate directions
543  for (unsigned i = 0; i < 2; i++)
544  {
545  mesh_velocity[i] += this->raw_dnodal_position_dt(l, i) * psif(l);
546  }
547  }
548  }
549 
550  // Get velocities and their derivatives from base flow problem
551  // -----------------------------------------------------------
552 
553  // Allocate storage for the velocity components of the base state
554  // solution (initialise to zero)
555  Vector<double> base_flow_u(3, 0.0);
556 
557  // Get the user-defined base state solution velocity components
558  get_base_flow_u(time, ipt, interpolated_x, base_flow_u);
559 
560  // Allocate storage for the derivatives of the base state solution's
561  // velocity components w.r.t. global coordinate (r and z)
562  // N.B. the derivatives of the base flow components w.r.t. the
563  // azimuthal coordinate direction (theta) are always zero since the
564  // base flow is axisymmetric
565  DenseMatrix<double> base_flow_dudx(3, 2, 0.0);
566 
567  // Get the derivatives of the user-defined base state solution
568  // velocity components w.r.t. global coordinates
569  get_base_flow_dudx(time, ipt, interpolated_x, base_flow_dudx);
570 
571  // Cache base flow velocities and their derivatives
572  const double interpolated_ur = base_flow_u[0];
573  const double interpolated_uz = base_flow_u[1];
574  const double interpolated_utheta = base_flow_u[2];
575  const double interpolated_durdr = base_flow_dudx(0, 0);
576  const double interpolated_durdz = base_flow_dudx(0, 1);
577  const double interpolated_duzdr = base_flow_dudx(1, 0);
578  const double interpolated_duzdz = base_flow_dudx(1, 1);
579  const double interpolated_duthetadr = base_flow_dudx(2, 0);
580  const double interpolated_duthetadz = base_flow_dudx(2, 1);
581 
582  // Cache r-component of position
583  const double r = interpolated_x[0];
584 
585  // Cache unknowns
586  const double interpolated_UC = interpolated_u[0];
587  const double interpolated_US = interpolated_u[1];
588  const double interpolated_WC = interpolated_u[2];
589  const double interpolated_WS = interpolated_u[3];
590  const double interpolated_VC = interpolated_u[4];
591  const double interpolated_VS = interpolated_u[5];
592  const double interpolated_PC = interpolated_p[0];
593  const double interpolated_PS = interpolated_p[1];
594 
595  // Cache derivatives of the unknowns
596  const double interpolated_dUCdr = interpolated_dudx(0, 0);
597  const double interpolated_dUCdz = interpolated_dudx(0, 1);
598  const double interpolated_dUSdr = interpolated_dudx(1, 0);
599  const double interpolated_dUSdz = interpolated_dudx(1, 1);
600  const double interpolated_dWCdr = interpolated_dudx(2, 0);
601  const double interpolated_dWCdz = interpolated_dudx(2, 1);
602  const double interpolated_dWSdr = interpolated_dudx(3, 0);
603  const double interpolated_dWSdz = interpolated_dudx(3, 1);
604  const double interpolated_dVCdr = interpolated_dudx(4, 0);
605  const double interpolated_dVCdz = interpolated_dudx(4, 1);
606  const double interpolated_dVSdr = interpolated_dudx(5, 0);
607  const double interpolated_dVSdz = interpolated_dudx(5, 1);
608 
609  // ==================
610  // MOMENTUM EQUATIONS
611  // ==================
612 
613  // Loop over the test functions
614  for (unsigned l = 0; l < n_node; l++)
615  {
616  // Cache test functions and their derivatives
617  const double testf_ = testf(l);
618  const double dtestfdr = dtestfdx(l, 0);
619  const double dtestfdz = dtestfdx(l, 1);
620 
621  // ---------------------------------------------
622  // FIRST (RADIAL) MOMENTUM EQUATION: COSINE PART
623  // ---------------------------------------------
624 
625  // Get local equation number of first velocity value at this node
626  local_eqn = nodal_local_eqn(l, u_nodal_index[0]);
627 
628  // If it's not a boundary condition
629  if (local_eqn >= 0)
630  {
631  // Pressure gradient term
632  residuals[local_eqn] += interpolated_PC * (testf_ + r * dtestfdr) * W;
633 
634  // Stress tensor terms
635  residuals[local_eqn] -= visc_ratio * r * (1.0 + Gamma[0]) *
636  interpolated_dUCdr * dtestfdr * W;
637 
638  residuals[local_eqn] -=
639  visc_ratio * r *
640  (interpolated_dUCdz + Gamma[0] * interpolated_dWCdr) * dtestfdz * W;
641 
642  residuals[local_eqn] +=
643  visc_ratio *
644  ((k * Gamma[0] * interpolated_dVSdr) -
645  (k * (2.0 + Gamma[0]) * interpolated_VS / r) -
646  ((1.0 + Gamma[0] + (k * k)) * interpolated_UC / r)) *
647  testf_ * W;
648 
649  // Inertial terms (du/dt)
650  residuals[local_eqn] -= scaled_re_st * r * dudt[0] * testf_ * W;
651 
652  // Inertial terms (convective)
653  residuals[local_eqn] -= scaled_re *
654  (r * interpolated_ur * interpolated_dUCdr +
655  r * interpolated_durdr * interpolated_UC +
656  k * interpolated_utheta * interpolated_US -
657  2 * interpolated_utheta * interpolated_VC +
658  r * interpolated_uz * interpolated_dUCdz +
659  r * interpolated_durdz * interpolated_WC) *
660  testf_ * W;
661 
662  // Mesh velocity terms
663  if (!ALE_is_disabled)
664  {
665  for (unsigned j = 0; j < 2; j++)
666  {
667  residuals[local_eqn] += scaled_re_st * r * mesh_velocity[j] *
668  interpolated_dudx(0, j) * testf_ * W;
669  }
670  }
671 
672  // Calculate the Jacobian
673  // ----------------------
674 
675  if (flag)
676  {
677  // Loop over the velocity shape functions again
678  for (unsigned l2 = 0; l2 < n_node; l2++)
679  {
680  // Cache velocity shape functions and their derivatives
681  const double psif_ = psif[l2];
682  const double dpsifdr = dpsifdx(l2, 0);
683  const double dpsifdz = dpsifdx(l2, 1);
684 
685  // Radial velocity component (cosine part) U_k^C
686  local_unknown = nodal_local_eqn(l2, u_nodal_index[0]);
687  if (local_unknown >= 0)
688  {
689  if (flag == 2)
690  {
691  // Add the mass matrix
692  mass_matrix(local_eqn, local_unknown) +=
693  scaled_re_st * r * psif_ * testf_ * W;
694  }
695 
696  // Add contributions to the Jacobian matrix
697 
698  // Stress tensor terms
699  jacobian(local_eqn, local_unknown) -=
700  visc_ratio * r * (1.0 + Gamma[0]) * dpsifdr * dtestfdr * W;
701 
702  jacobian(local_eqn, local_unknown) -=
703  visc_ratio * r * dpsifdz * dtestfdz * W;
704 
705  jacobian(local_eqn, local_unknown) -=
706  visc_ratio * (1.0 + Gamma[0] + (k * k)) * psif_ * testf_ * W /
707  r;
708 
709  // Inertial terms (du/dt)
710  jacobian(local_eqn, local_unknown) -=
711  scaled_re_st * r *
712  node_pt(l2)->time_stepper_pt()->weight(1, 0) * psif_ *
713  testf_ * W;
714 
715  // Inertial terms (convective)
716  jacobian(local_eqn, local_unknown) -=
717  scaled_re * r *
718  (psif_ * interpolated_durdr + interpolated_ur * dpsifdr +
719  interpolated_uz * dpsifdz) *
720  testf_ * W;
721 
722  // Mesh velocity terms
723  if (!ALE_is_disabled)
724  {
725  for (unsigned j = 0; j < 2; j++)
726  {
727  jacobian(local_eqn, local_unknown) +=
728  scaled_re_st * r * mesh_velocity[j] * dpsifdx(l2, j) *
729  testf_ * W;
730  }
731  }
732  }
733 
734  // Radial velocity component (sine part) U_k^S
735  local_unknown = nodal_local_eqn(l2, u_nodal_index[1]);
736  if (local_unknown >= 0)
737  {
738  // Convective term
739  jacobian(local_eqn, local_unknown) -=
740  scaled_re * k * interpolated_utheta * psif_ * testf_ * W;
741  }
742 
743  // Axial velocity component (cosine part) W_k^C
744  local_unknown = nodal_local_eqn(l2, u_nodal_index[2]);
745  if (local_unknown >= 0)
746  {
747  // Stress tensor term
748  jacobian(local_eqn, local_unknown) -=
749  visc_ratio * Gamma[0] * r * dpsifdr * dtestfdz * W;
750 
751  // Convective term
752  jacobian(local_eqn, local_unknown) -=
753  scaled_re * r * interpolated_durdz * psif_ * testf_ * W;
754  }
755 
756  // Axial velocity component (sine part) W_k^S
757  // has no contribution
758 
759  // Azimuthal velocity component (cosine part) V_k^C
760  local_unknown = nodal_local_eqn(l2, u_nodal_index[4]);
761  if (local_unknown >= 0)
762  {
763  // Convective term
764  jacobian(local_eqn, local_unknown) +=
765  scaled_re * 2 * interpolated_utheta * psif_ * testf_ * W;
766  }
767 
768  // Azimuthal velocity component (sine part) V_k^S
769  local_unknown = nodal_local_eqn(l2, u_nodal_index[5]);
770  if (local_unknown >= 0)
771  {
772  // Stress tensor term
773  jacobian(local_eqn, local_unknown) +=
774  visc_ratio *
775  ((Gamma[0] * k * dpsifdr) -
776  (k * (2.0 + Gamma[0]) * psif_ / r)) *
777  testf_ * W;
778  }
779  } // End of loop over velocity shape functions
780 
781  // Now loop over pressure shape functions
782  // (This is the contribution from pressure gradient)
783  for (unsigned l2 = 0; l2 < n_pres; l2++)
784  {
785  // Cosine part P_k^C
786  local_unknown = p_local_eqn(l2, 0);
787  if (local_unknown >= 0)
788  {
789  jacobian(local_eqn, local_unknown) +=
790  psip[l2] * (testf_ + r * dtestfdr) * W;
791  }
792 
793  // Sine part P_k^S has no contribution
794 
795  } // End of loop over pressure shape functions
796  } // End of Jacobian calculation
797 
798  } // End of if not boundary condition statement
799 
800  // --------------------------------------------
801  // SECOND (RADIAL) MOMENTUM EQUATION: SINE PART
802  // --------------------------------------------
803 
804  // Get local equation number of second velocity value at this node
805  local_eqn = nodal_local_eqn(l, u_nodal_index[1]);
806 
807  // If it's not a boundary condition
808  if (local_eqn >= 0)
809  {
810  // Pressure gradient term
811  residuals[local_eqn] += interpolated_PS * (testf_ + r * dtestfdr) * W;
812 
813  // Stress tensor terms
814  residuals[local_eqn] -= visc_ratio * r * (1.0 + Gamma[0]) *
815  interpolated_dUSdr * dtestfdr * W;
816 
817  residuals[local_eqn] -=
818  visc_ratio * r *
819  (interpolated_dUSdz + Gamma[0] * interpolated_dWSdr) * dtestfdz * W;
820 
821  residuals[local_eqn] -=
822  visc_ratio *
823  ((k * Gamma[0] * interpolated_dVCdr) -
824  (k * (2.0 + Gamma[0]) * interpolated_VC / r) +
825  ((1.0 + Gamma[0] + (k * k)) * interpolated_US / r)) *
826  testf_ * W;
827 
828  // Inertial terms (du/dt)
829  residuals[local_eqn] -= scaled_re_st * r * dudt[1] * testf_ * W;
830 
831  // Inertial terms (convective)
832  residuals[local_eqn] -= scaled_re *
833  (r * interpolated_ur * interpolated_dUSdr +
834  r * interpolated_durdr * interpolated_US -
835  k * interpolated_utheta * interpolated_UC -
836  2 * interpolated_utheta * interpolated_VS +
837  r * interpolated_uz * interpolated_dUSdz +
838  r * interpolated_durdz * interpolated_WS) *
839  testf_ * W;
840 
841  // Mesh velocity terms
842  if (!ALE_is_disabled)
843  {
844  for (unsigned j = 0; j < 2; j++)
845  {
846  residuals[local_eqn] += scaled_re_st * r * mesh_velocity[j] *
847  interpolated_dudx(1, j) * testf_ * W;
848  }
849  }
850 
851  // Calculate the Jacobian
852  // ----------------------
853 
854  if (flag)
855  {
856  // Loop over the velocity shape functions again
857  for (unsigned l2 = 0; l2 < n_node; l2++)
858  {
859  // Cache velocity shape functions and their derivatives
860  const double psif_ = psif[l2];
861  const double dpsifdr = dpsifdx(l2, 0);
862  const double dpsifdz = dpsifdx(l2, 1);
863 
864  // Radial velocity component (cosine part) U_k^C
865  local_unknown = nodal_local_eqn(l2, u_nodal_index[0]);
866  if (local_unknown >= 0)
867  {
868  // Convective term
869  jacobian(local_eqn, local_unknown) +=
870  scaled_re * k * interpolated_utheta * psif_ * testf_ * W;
871  }
872 
873  // Radial velocity component (sine part) U_k^S
874  local_unknown = nodal_local_eqn(l2, u_nodal_index[1]);
875  if (local_unknown >= 0)
876  {
877  if (flag == 2)
878  {
879  // Add the mass matrix
880  mass_matrix(local_eqn, local_unknown) +=
881  scaled_re_st * r * psif_ * testf_ * W;
882  }
883 
884  // Add contributions to the Jacobian matrix
885 
886  // Stress tensor terms
887  jacobian(local_eqn, local_unknown) -=
888  visc_ratio * r * (1.0 + Gamma[0]) * dpsifdr * dtestfdr * W;
889 
890  jacobian(local_eqn, local_unknown) -=
891  visc_ratio * r * dpsifdz * dtestfdz * W;
892 
893  jacobian(local_eqn, local_unknown) -=
894  visc_ratio * (1.0 + Gamma[0] + (k * k)) * psif_ * testf_ * W /
895  r;
896 
897  // Inertial terms (du/dt)
898  jacobian(local_eqn, local_unknown) -=
899  scaled_re_st * r *
900  node_pt(l2)->time_stepper_pt()->weight(1, 0) * psif_ *
901  testf_ * W;
902 
903  // Inertial terms (convective)
904  jacobian(local_eqn, local_unknown) -=
905  scaled_re * r *
906  (psif_ * interpolated_durdr + interpolated_ur * dpsifdr +
907  interpolated_uz * dpsifdz) *
908  testf_ * W;
909 
910  // Mesh velocity terms
911  if (!ALE_is_disabled)
912  {
913  for (unsigned j = 0; j < 2; j++)
914  {
915  jacobian(local_eqn, local_unknown) +=
916  scaled_re_st * r * mesh_velocity[j] * dpsifdx(l2, j) *
917  testf_ * W;
918  }
919  }
920  }
921 
922  // Axial velocity component (cosine part) W_k^C
923  // has no contribution
924 
925  // Axial velocity component (sine part) W_k^S
926  local_unknown = nodal_local_eqn(l2, u_nodal_index[3]);
927  if (local_unknown >= 0)
928  {
929  // Stress tensor term
930  jacobian(local_eqn, local_unknown) -=
931  visc_ratio * Gamma[0] * r * dpsifdr * dtestfdz * W;
932 
933  // Convective term
934  jacobian(local_eqn, local_unknown) -=
935  scaled_re * r * interpolated_durdz * psif_ * testf_ * W;
936  }
937 
938  // Azimuthal velocity component (cosine part) V_k^C
939  local_unknown = nodal_local_eqn(l2, u_nodal_index[4]);
940  if (local_unknown >= 0)
941  {
942  // Stress tensor terms
943  jacobian(local_eqn, local_unknown) -=
944  visc_ratio *
945  ((Gamma[0] * k * dpsifdr) -
946  (k * (2.0 + Gamma[0]) * psif_ / r)) *
947  testf_ * W;
948  }
949 
950  // Azimuthal velocity component (sine part) V_k^S
951  local_unknown = nodal_local_eqn(l2, u_nodal_index[5]);
952  if (local_unknown >= 0)
953  {
954  // Convective term
955  jacobian(local_eqn, local_unknown) +=
956  scaled_re * 2 * interpolated_utheta * psif_ * testf_ * W;
957  }
958  } // End of loop over velocity shape functions
959 
960  // Now loop over pressure shape functions
961  // (This is the contribution from pressure gradient)
962  for (unsigned l2 = 0; l2 < n_pres; l2++)
963  {
964  // Cosine part P_k^C has no contribution
965 
966  // Sine part P_k^S
967  local_unknown = p_local_eqn(l2, 1);
968  if (local_unknown >= 0)
969  {
970  jacobian(local_eqn, local_unknown) +=
971  psip[l2] * (testf_ + r * dtestfdr) * W;
972  }
973  } // End of loop over pressure shape functions
974  } // End of Jacobian calculation
975 
976  } // End of if not boundary condition statement
977 
978  // --------------------------------------------
979  // THIRD (AXIAL) MOMENTUM EQUATION: COSINE PART
980  // --------------------------------------------
981 
982  // Get local equation number of third velocity value at this node
983  local_eqn = nodal_local_eqn(l, u_nodal_index[2]);
984 
985  // If it's not a boundary condition
986  if (local_eqn >= 0)
987  {
988  // Pressure gradient term
989  residuals[local_eqn] += r * interpolated_PC * dtestfdz * W;
990 
991  // Stress tensor terms
992  residuals[local_eqn] -=
993  visc_ratio * r *
994  (interpolated_dWCdr + Gamma[1] * interpolated_dUCdz) * dtestfdr * W;
995 
996  residuals[local_eqn] -= visc_ratio * r * (1.0 + Gamma[1]) *
997  interpolated_dWCdz * dtestfdz * W;
998 
999  residuals[local_eqn] +=
1000  visc_ratio * k *
1001  (Gamma[1] * interpolated_dVSdz - k * interpolated_WC / r) * testf_ *
1002  W;
1003 
1004  // Inertial terms (du/dt)
1005  residuals[local_eqn] -= scaled_re_st * r * dudt[2] * testf_ * W;
1006 
1007  // Inertial terms (convective)
1008  residuals[local_eqn] -= scaled_re *
1009  (r * interpolated_ur * interpolated_dWCdr +
1010  r * interpolated_duzdr * interpolated_UC +
1011  k * interpolated_utheta * interpolated_WS +
1012  r * interpolated_uz * interpolated_dWCdz +
1013  r * interpolated_duzdz * interpolated_WC) *
1014  testf_ * W;
1015 
1016  // Mesh velocity terms
1017  if (!ALE_is_disabled)
1018  {
1019  for (unsigned j = 0; j < 2; j++)
1020  {
1021  residuals[local_eqn] += scaled_re_st * r * mesh_velocity[j] *
1022  interpolated_dudx(2, j) * testf_ * W;
1023  }
1024  }
1025 
1026  // Calculate the Jacobian
1027  // ----------------------
1028 
1029  if (flag)
1030  {
1031  // Loop over the velocity shape functions again
1032  for (unsigned l2 = 0; l2 < n_node; l2++)
1033  {
1034  // Cache velocity shape functions and their derivatives
1035  const double psif_ = psif[l2];
1036  const double dpsifdr = dpsifdx(l2, 0);
1037  const double dpsifdz = dpsifdx(l2, 1);
1038 
1039  // Radial velocity component (cosine part) U_k^C
1040  local_unknown = nodal_local_eqn(l2, u_nodal_index[0]);
1041  if (local_unknown >= 0)
1042  {
1043  // Stress tensor term
1044  jacobian(local_eqn, local_unknown) -=
1045  visc_ratio * r * Gamma[1] * dpsifdz * dtestfdr * W;
1046 
1047  // Convective term
1048  jacobian(local_eqn, local_unknown) -=
1049  scaled_re * r * psif_ * interpolated_duzdr * testf_ * W;
1050  }
1051 
1052  // Radial velocity component (sine part) U_k^S
1053  // has no contribution
1054 
1055  // Axial velocity component (cosine part) W_k^C
1056  local_unknown = nodal_local_eqn(l2, u_nodal_index[2]);
1057  if (local_unknown >= 0)
1058  {
1059  if (flag == 2)
1060  {
1061  // Add the mass matrix
1062  mass_matrix(local_eqn, local_unknown) +=
1063  scaled_re_st * r * psif_ * testf_ * W;
1064  }
1065 
1066  // Add contributions to the Jacobian matrix
1067 
1068  // Stress tensor terms
1069  jacobian(local_eqn, local_unknown) -=
1070  visc_ratio * r * dpsifdr * dtestfdr * W;
1071 
1072  jacobian(local_eqn, local_unknown) -=
1073  visc_ratio * r * (1.0 + Gamma[1]) * dpsifdz * dtestfdz * W;
1074 
1075  jacobian(local_eqn, local_unknown) -=
1076  visc_ratio * k * k * psif_ * testf_ * W / r;
1077 
1078  // Inertial terms (du/dt)
1079  jacobian(local_eqn, local_unknown) -=
1080  scaled_re_st * r *
1081  node_pt(l2)->time_stepper_pt()->weight(1, 0) * psif_ *
1082  testf_ * W;
1083 
1084  // Inertial terms (convective)
1085  jacobian(local_eqn, local_unknown) -=
1086  scaled_re * r *
1087  (interpolated_ur * dpsifdr + psif_ * interpolated_duzdz +
1088  interpolated_uz * dpsifdz) *
1089  testf_ * W;
1090 
1091 
1092  // Mesh velocity terms
1093  if (!ALE_is_disabled)
1094  {
1095  for (unsigned j = 0; j < 2; j++)
1096  {
1097  jacobian(local_eqn, local_unknown) +=
1098  scaled_re_st * r * mesh_velocity[j] * dpsifdx(l2, j) *
1099  testf_ * W;
1100  }
1101  }
1102  }
1103 
1104  // Axial velocity component (sine part) W_k^S
1105  local_unknown = nodal_local_eqn(l2, u_nodal_index[3]);
1106  if (local_unknown >= 0)
1107  {
1108  // Convective term
1109  jacobian(local_eqn, local_unknown) -=
1110  scaled_re * k * interpolated_utheta * psif_ * testf_ * W;
1111  }
1112 
1113  // Azimuthal velocity component (cosine part) V_k^C
1114  // has no contribution
1115 
1116  // Azimuthal velocity component (sine part) V_k^S
1117  local_unknown = nodal_local_eqn(l2, u_nodal_index[5]);
1118  if (local_unknown >= 0)
1119  {
1120  // Stress tensor term
1121  jacobian(local_eqn, local_unknown) +=
1122  visc_ratio * Gamma[1] * k * dpsifdz * testf_ * W;
1123  }
1124  } // End of loop over velocity shape functions
1125 
1126  // Now loop over pressure shape functions
1127  // (This is the contribution from pressure gradient)
1128  for (unsigned l2 = 0; l2 < n_pres; l2++)
1129  {
1130  // Cosine part P_k^C
1131  local_unknown = p_local_eqn(l2, 0);
1132  if (local_unknown >= 0)
1133  {
1134  jacobian(local_eqn, local_unknown) +=
1135  r * psip[l2] * dtestfdz * W;
1136  }
1137 
1138  // Sine part P_k^S has no contribution
1139 
1140  } // End of loop over pressure shape functions
1141  } // End of Jacobian calculation
1142 
1143  } // End of if not boundary condition statement
1144 
1145  // -------------------------------------------
1146  // FOURTH (AXIAL) MOMENTUM EQUATION: SINE PART
1147  // -------------------------------------------
1148 
1149  // Get local equation number of fourth velocity value at this node
1150  local_eqn = nodal_local_eqn(l, u_nodal_index[3]);
1151 
1152  // If it's not a boundary condition
1153  if (local_eqn >= 0)
1154  {
1155  // Pressure gradient term
1156  residuals[local_eqn] += r * interpolated_PS * dtestfdz * W;
1157 
1158  // Stress tensor terms
1159  residuals[local_eqn] -=
1160  visc_ratio * r *
1161  (interpolated_dWSdr + Gamma[1] * interpolated_dUSdz) * dtestfdr * W;
1162 
1163  residuals[local_eqn] -= visc_ratio * r * (1.0 + Gamma[1]) *
1164  interpolated_dWSdz * dtestfdz * W;
1165 
1166  residuals[local_eqn] -=
1167  visc_ratio * k *
1168  (Gamma[1] * interpolated_dVCdz + k * interpolated_WS / r) * testf_ *
1169  W;
1170 
1171  // Inertial terms (du/dt)
1172  residuals[local_eqn] -= scaled_re_st * r * dudt[3] * testf_ * W;
1173 
1174  // Inertial terms (convective)
1175  residuals[local_eqn] -= scaled_re *
1176  (r * interpolated_ur * interpolated_dWSdr +
1177  r * interpolated_duzdr * interpolated_US -
1178  k * interpolated_utheta * interpolated_WC +
1179  r * interpolated_uz * interpolated_dWSdz +
1180  r * interpolated_duzdz * interpolated_WS) *
1181  testf_ * W;
1182 
1183  // Mesh velocity terms
1184  if (!ALE_is_disabled)
1185  {
1186  for (unsigned j = 0; j < 2; j++)
1187  {
1188  residuals[local_eqn] += scaled_re_st * r * mesh_velocity[j] *
1189  interpolated_dudx(3, j) * testf_ * W;
1190  }
1191  }
1192 
1193  // Calculate the Jacobian
1194  // ----------------------
1195 
1196  if (flag)
1197  {
1198  // Loop over the velocity shape functions again
1199  for (unsigned l2 = 0; l2 < n_node; l2++)
1200  {
1201  // Cache velocity shape functions and their derivatives
1202  const double psif_ = psif[l2];
1203  const double dpsifdr = dpsifdx(l2, 0);
1204  const double dpsifdz = dpsifdx(l2, 1);
1205 
1206  // Radial velocity component (cosine part) U_k^C
1207  // has no contribution
1208 
1209  // Radial velocity component (sine part) U_k^S
1210  local_unknown = nodal_local_eqn(l2, u_nodal_index[1]);
1211  if (local_unknown >= 0)
1212  {
1213  // Stress tensor term
1214  jacobian(local_eqn, local_unknown) -=
1215  visc_ratio * r * Gamma[1] * dpsifdz * dtestfdr * W;
1216 
1217  // Convective term
1218  jacobian(local_eqn, local_unknown) -=
1219  scaled_re * r * psif_ * interpolated_duzdr * testf_ * W;
1220  }
1221 
1222  // Axial velocity component (cosine part) W_k^S
1223  local_unknown = nodal_local_eqn(l2, u_nodal_index[2]);
1224  if (local_unknown >= 0)
1225  {
1226  // Convective term
1227  jacobian(local_eqn, local_unknown) +=
1228  scaled_re * k * interpolated_utheta * psif_ * testf_ * W;
1229  }
1230 
1231  // Axial velocity component (sine part) W_k^S
1232  local_unknown = nodal_local_eqn(l2, u_nodal_index[3]);
1233  if (local_unknown >= 0)
1234  {
1235  if (flag == 2)
1236  {
1237  // Add the mass matrix
1238  mass_matrix(local_eqn, local_unknown) +=
1239  scaled_re_st * r * psif_ * testf_ * W;
1240  }
1241 
1242  // Add contributions to the Jacobian matrix
1243 
1244  // Stress tensor terms
1245  jacobian(local_eqn, local_unknown) -=
1246  visc_ratio * r * dpsifdr * dtestfdr * W;
1247 
1248  jacobian(local_eqn, local_unknown) -=
1249  visc_ratio * r * (1.0 + Gamma[1]) * dpsifdz * dtestfdz * W;
1250 
1251  jacobian(local_eqn, local_unknown) -=
1252  visc_ratio * k * k * psif_ * testf_ * W / r;
1253 
1254  // Inertial terms (du/dt)
1255  jacobian(local_eqn, local_unknown) -=
1256  scaled_re_st * r *
1257  node_pt(l2)->time_stepper_pt()->weight(1, 0) * psif_ *
1258  testf_ * W;
1259 
1260  // Inertial terms (convective)
1261  jacobian(local_eqn, local_unknown) -=
1262  scaled_re * r *
1263  (interpolated_ur * dpsifdr + psif_ * interpolated_duzdz +
1264  interpolated_uz * dpsifdz) *
1265  testf_ * W;
1266 
1267 
1268  // Mesh velocity terms
1269  if (!ALE_is_disabled)
1270  {
1271  for (unsigned j = 0; j < 2; j++)
1272  {
1273  jacobian(local_eqn, local_unknown) +=
1274  scaled_re_st * r * mesh_velocity[j] * dpsifdx(l2, j) *
1275  testf_ * W;
1276  }
1277  }
1278  }
1279 
1280  // Azimuthal velocity component (cosine part) V_k^C
1281  local_unknown = nodal_local_eqn(l2, u_nodal_index[4]);
1282  if (local_unknown >= 0)
1283  {
1284  // Stress tensor term
1285  jacobian(local_eqn, local_unknown) -=
1286  visc_ratio * Gamma[1] * k * dpsifdz * testf_ * W;
1287  }
1288 
1289  // Azimuthal velocity component (sine part) V_k^S
1290  // has no contribution
1291 
1292  } // End of loop over velocity shape functions
1293 
1294  // Now loop over pressure shape functions
1295  // (This is the contribution from pressure gradient)
1296  for (unsigned l2 = 0; l2 < n_pres; l2++)
1297  {
1298  // Cosine part P_k^C has no contribution
1299 
1300  // Sine part P_k^S
1301  local_unknown = p_local_eqn(l2, 1);
1302  if (local_unknown >= 0)
1303  {
1304  jacobian(local_eqn, local_unknown) +=
1305  r * psip[l2] * dtestfdz * W;
1306  }
1307  } // End of loop over pressure shape functions
1308  } // End of Jacobian calculation
1309 
1310  } // End of if not boundary condition statement
1311 
1312  // ------------------------------------------------
1313  // FIFTH (AZIMUTHAL) MOMENTUM EQUATION: COSINE PART
1314  // ------------------------------------------------
1315 
1316  // Get local equation number of fifth velocity value at this node
1317  local_eqn = nodal_local_eqn(l, u_nodal_index[4]);
1318 
1319  // If it's not a boundary condition
1320  if (local_eqn >= 0)
1321  {
1322  // Pressure gradient term
1323  residuals[local_eqn] -= k * interpolated_PS * testf_ * W;
1324 
1325  // Stress tensor terms
1326  residuals[local_eqn] +=
1327  visc_ratio *
1328  (-r * interpolated_dVCdr - k * Gamma[0] * interpolated_US +
1329  Gamma[0] * interpolated_VC) *
1330  dtestfdr * W;
1331 
1332  residuals[local_eqn] -=
1333  visc_ratio *
1334  (k * Gamma[0] * interpolated_WS + r * interpolated_dVCdz) *
1335  dtestfdz * W;
1336 
1337  residuals[local_eqn] +=
1338  visc_ratio *
1339  (Gamma[0] * interpolated_dVCdr +
1340  k * (2.0 + Gamma[0]) * interpolated_US / r -
1341  (1.0 + (k * k) + (k * k * Gamma[0])) * interpolated_VC / r) *
1342  testf_ * W;
1343 
1344  // Inertial terms (du/dt)
1345  residuals[local_eqn] -= scaled_re_st * r * dudt[4] * testf_ * W;
1346 
1347  // Inertial terms (convective)
1348  residuals[local_eqn] -=
1349  scaled_re *
1350  (r * interpolated_ur * interpolated_dVCdr +
1351  r * interpolated_duthetadr * interpolated_UC +
1352  k * interpolated_utheta * interpolated_VS +
1353  interpolated_utheta * interpolated_UC +
1354  interpolated_ur * interpolated_VC +
1355  r * interpolated_uz * interpolated_dVCdz +
1356  r * interpolated_duthetadz * interpolated_WC) *
1357  testf_ * W;
1358 
1359  // Mesh velocity terms
1360  if (!ALE_is_disabled)
1361  {
1362  for (unsigned j = 0; j < 2; j++)
1363  {
1364  residuals[local_eqn] += scaled_re_st * r * mesh_velocity[j] *
1365  interpolated_dudx(4, j) * testf_ * W;
1366  }
1367  }
1368 
1369  // Calculate the Jacobian
1370  // ----------------------
1371 
1372  if (flag)
1373  {
1374  // Loop over the velocity shape functions again
1375  for (unsigned l2 = 0; l2 < n_node; l2++)
1376  {
1377  // Cache velocity shape functions and their derivatives
1378  const double psif_ = psif[l2];
1379  const double dpsifdr = dpsifdx(l2, 0);
1380  const double dpsifdz = dpsifdx(l2, 1);
1381 
1382  // Radial velocity component (cosine part) U_k^C
1383  local_unknown = nodal_local_eqn(l2, u_nodal_index[0]);
1384  if (local_unknown >= 0)
1385  {
1386  // Convective terms
1387  jacobian(local_eqn, local_unknown) -=
1388  scaled_re *
1389  (r * interpolated_duthetadr + interpolated_utheta) * psif_ *
1390  testf_ * W;
1391  }
1392 
1393  // Radial velocity component (sine part) U_k^S
1394  local_unknown = nodal_local_eqn(l2, u_nodal_index[1]);
1395  if (local_unknown >= 0)
1396  {
1397  // Stress tensor terms
1398  jacobian(local_eqn, local_unknown) +=
1399  visc_ratio * k * psif_ *
1400  (((2.0 + Gamma[0]) * testf_ / r) - (Gamma[0] * dtestfdr)) * W;
1401  }
1402 
1403  // Axial velocity component (cosine part) W_k^C
1404  local_unknown = nodal_local_eqn(l2, u_nodal_index[2]);
1405  if (local_unknown >= 0)
1406  {
1407  // Convective term
1408  jacobian(local_eqn, local_unknown) -=
1409  scaled_re * r * psif_ * interpolated_duthetadz * testf_ * W;
1410  }
1411 
1412  // Axial velocity component (sine part) W_k^S
1413  local_unknown = nodal_local_eqn(l2, u_nodal_index[3]);
1414  if (local_unknown >= 0)
1415  {
1416  // Stress tensor term
1417  jacobian(local_eqn, local_unknown) -=
1418  visc_ratio * k * Gamma[0] * psif_ * dtestfdz * W;
1419  }
1420 
1421  // Azimuthal velocity component (cosine part) V_k^C
1422  local_unknown = nodal_local_eqn(l2, u_nodal_index[4]);
1423  if (local_unknown >= 0)
1424  {
1425  if (flag == 2)
1426  {
1427  // Add the mass matrix
1428  mass_matrix(local_eqn, local_unknown) +=
1429  scaled_re_st * r * psif_ * testf_ * W;
1430  }
1431 
1432  // Add contributions to the Jacobian matrix
1433 
1434  // Stress tensor terms
1435  jacobian(local_eqn, local_unknown) +=
1436  visc_ratio * (Gamma[0] * psif_ - r * dpsifdr) * dtestfdr * W;
1437 
1438  jacobian(local_eqn, local_unknown) -=
1439  visc_ratio * r * dpsifdz * dtestfdz * W;
1440 
1441  jacobian(local_eqn, local_unknown) +=
1442  visc_ratio *
1443  (Gamma[0] * dpsifdr -
1444  (1.0 + (k * k) + (k * k * Gamma[0])) * psif_ / r) *
1445  testf_ * W;
1446 
1447  // Inertial terms (du/dt)
1448  jacobian(local_eqn, local_unknown) -=
1449  scaled_re_st * r *
1450  node_pt(l2)->time_stepper_pt()->weight(1, 0) * psif_ *
1451  testf_ * W;
1452 
1453  // Inertial terms (convective)
1454  jacobian(local_eqn, local_unknown) -=
1455  scaled_re *
1456  (r * interpolated_ur * dpsifdr + interpolated_ur * psif_ +
1457  r * interpolated_uz * dpsifdz) *
1458  testf_ * W;
1459 
1460  // Mesh velocity terms
1461  if (!ALE_is_disabled)
1462  {
1463  for (unsigned j = 0; j < 2; j++)
1464  {
1465  jacobian(local_eqn, local_unknown) +=
1466  scaled_re_st * r * mesh_velocity[j] * dpsifdx(l2, j) *
1467  testf_ * W;
1468  }
1469  }
1470  }
1471 
1472  // Azimuthal velocity component (sine part) V_k^S
1473  local_unknown = nodal_local_eqn(l2, u_nodal_index[5]);
1474  if (local_unknown >= 0)
1475  {
1476  // Convective term
1477  jacobian(local_eqn, local_unknown) -=
1478  scaled_re * k * interpolated_utheta * psif_ * testf_ * W;
1479  }
1480  } // End of loop over velocity shape functions
1481 
1482  // Now loop over pressure shape functions
1483  // (This is the contribution from pressure gradient)
1484  for (unsigned l2 = 0; l2 < n_pres; l2++)
1485  {
1486  // Cosine part P_k^C has no contribution
1487 
1488  // Sine part P_k^S
1489  local_unknown = p_local_eqn(l2, 1);
1490  if (local_unknown >= 0)
1491  {
1492  jacobian(local_eqn, local_unknown) -= k * psip[l2] * testf_ * W;
1493  }
1494  } // End of loop over pressure shape functions
1495  } // End of Jacobian calculation
1496 
1497  } // End of if not boundary condition statement
1498 
1499  // ----------------------------------------------
1500  // SIXTH (AZIMUTHAL) MOMENTUM EQUATION: SINE PART
1501  // ----------------------------------------------
1502 
1503  // Get local equation number of sixth velocity value at this node
1504  local_eqn = nodal_local_eqn(l, u_nodal_index[5]);
1505 
1506  // If it's not a boundary condition
1507  if (local_eqn >= 0)
1508  {
1509  // Pressure gradient term
1510  residuals[local_eqn] += k * interpolated_PC * testf_ * W;
1511 
1512  // Stress tensor terms
1513  residuals[local_eqn] +=
1514  visc_ratio *
1515  (-r * interpolated_dVSdr + k * Gamma[0] * interpolated_UC +
1516  Gamma[0] * interpolated_VS) *
1517  dtestfdr * W;
1518 
1519  residuals[local_eqn] +=
1520  visc_ratio *
1521  (k * Gamma[0] * interpolated_WC - r * interpolated_dVSdz) *
1522  dtestfdz * W;
1523 
1524  residuals[local_eqn] +=
1525  visc_ratio *
1526  (Gamma[0] * interpolated_dVSdr -
1527  k * (2.0 + Gamma[0]) * interpolated_UC / r -
1528  (1.0 + (k * k) + (k * k * Gamma[0])) * interpolated_VS / r) *
1529  testf_ * W;
1530 
1531  // Inertial terms (du/dt)
1532  residuals[local_eqn] -= scaled_re_st * r * dudt[5] * testf_ * W;
1533 
1534  // Inertial terms (convective)
1535  residuals[local_eqn] -=
1536  scaled_re *
1537  (r * interpolated_ur * interpolated_dVSdr +
1538  r * interpolated_duthetadr * interpolated_US -
1539  k * interpolated_utheta * interpolated_VC +
1540  interpolated_utheta * interpolated_US +
1541  interpolated_ur * interpolated_VS +
1542  r * interpolated_uz * interpolated_dVSdz +
1543  r * interpolated_duthetadz * interpolated_WS) *
1544  testf_ * W;
1545 
1546  // Mesh velocity terms
1547  if (!ALE_is_disabled)
1548  {
1549  for (unsigned j = 0; j < 2; j++)
1550  {
1551  residuals[local_eqn] += scaled_re_st * r * mesh_velocity[j] *
1552  interpolated_dudx(5, j) * testf_ * W;
1553  }
1554  }
1555 
1556  // Calculate the Jacobian
1557  // ----------------------
1558 
1559  if (flag)
1560  {
1561  // Loop over the velocity shape functions again
1562  for (unsigned l2 = 0; l2 < n_node; l2++)
1563  {
1564  // Cache velocity shape functions and their derivatives
1565  const double psif_ = psif[l2];
1566  const double dpsifdr = dpsifdx(l2, 0);
1567  const double dpsifdz = dpsifdx(l2, 1);
1568 
1569  // Radial velocity component (cosine part) U_k^C
1570  local_unknown = nodal_local_eqn(l2, u_nodal_index[0]);
1571  if (local_unknown >= 0)
1572  {
1573  // Stress tensor terms
1574  jacobian(local_eqn, local_unknown) +=
1575  visc_ratio * k * psif_ *
1576  ((Gamma[0] * dtestfdr) - ((2.0 + Gamma[0]) * testf_ / r)) * W;
1577  }
1578 
1579  // Radial velocity component (sine part) U_k^S
1580  local_unknown = nodal_local_eqn(l2, u_nodal_index[1]);
1581  if (local_unknown >= 0)
1582  {
1583  // Convective terms
1584  jacobian(local_eqn, local_unknown) -=
1585  scaled_re *
1586  (r * interpolated_duthetadr + interpolated_utheta) * psif_ *
1587  testf_ * W;
1588  }
1589 
1590  // Axial velocity component (cosine part) W_k^C
1591  local_unknown = nodal_local_eqn(l2, u_nodal_index[2]);
1592  if (local_unknown >= 0)
1593  {
1594  // Stress tensor term
1595  jacobian(local_eqn, local_unknown) +=
1596  visc_ratio * k * Gamma[0] * psif_ * dtestfdz * W;
1597  }
1598 
1599  // Axial velocity component (sine part) W_k^S
1600  local_unknown = nodal_local_eqn(l2, u_nodal_index[3]);
1601  if (local_unknown >= 0)
1602  {
1603  // Convective term
1604  jacobian(local_eqn, local_unknown) -=
1605  scaled_re * r * psif_ * interpolated_duthetadz * testf_ * W;
1606  }
1607 
1608  // Azimuthal velocity component (cosine part) V_k^C
1609  local_unknown = nodal_local_eqn(l2, u_nodal_index[4]);
1610  if (local_unknown >= 0)
1611  {
1612  // Convective term
1613  jacobian(local_eqn, local_unknown) +=
1614  scaled_re * k * interpolated_utheta * psif_ * testf_ * W;
1615  }
1616 
1617  // Azimuthal velocity component (sine part) V_k^S
1618  local_unknown = nodal_local_eqn(l2, u_nodal_index[5]);
1619  if (local_unknown >= 0)
1620  {
1621  if (flag == 2)
1622  {
1623  // Add the mass matrix
1624  mass_matrix(local_eqn, local_unknown) +=
1625  scaled_re_st * r * psif_ * testf_ * W;
1626  }
1627 
1628  // Add contributions to the Jacobian matrix
1629 
1630  // Stress tensor terms
1631  jacobian(local_eqn, local_unknown) +=
1632  visc_ratio * (Gamma[0] * psif_ - r * dpsifdr) * dtestfdr * W;
1633 
1634  jacobian(local_eqn, local_unknown) -=
1635  visc_ratio * r * dpsifdz * dtestfdz * W;
1636 
1637  jacobian(local_eqn, local_unknown) +=
1638  visc_ratio *
1639  (Gamma[0] * dpsifdr -
1640  (1.0 + (k * k) + (k * k * Gamma[0])) * psif_ / r) *
1641  testf_ * W;
1642 
1643  // Inertial terms (du/dt)
1644  jacobian(local_eqn, local_unknown) -=
1645  scaled_re_st * r *
1646  node_pt(l2)->time_stepper_pt()->weight(1, 0) * psif_ *
1647  testf_ * W;
1648 
1649  // Convective terms
1650  jacobian(local_eqn, local_unknown) -=
1651  scaled_re *
1652  (r * interpolated_ur * dpsifdr + interpolated_ur * psif_ +
1653  r * interpolated_uz * dpsifdz) *
1654  testf_ * W;
1655 
1656  // Mesh velocity terms
1657  if (!ALE_is_disabled)
1658  {
1659  for (unsigned j = 0; j < 2; j++)
1660  {
1661  jacobian(local_eqn, local_unknown) +=
1662  scaled_re_st * r * mesh_velocity[j] * dpsifdx(l2, j) *
1663  testf_ * W;
1664  }
1665  }
1666  }
1667  } // End of loop over velocity shape functions
1668 
1669  // Now loop over pressure shape functions
1670  // (This is the contribution from pressure gradient)
1671  for (unsigned l2 = 0; l2 < n_pres; l2++)
1672  {
1673  // Cosine part P_k^C
1674  local_unknown = p_local_eqn(l2, 0);
1675  if (local_unknown >= 0)
1676  {
1677  jacobian(local_eqn, local_unknown) += k * psip[l2] * testf_ * W;
1678  }
1679 
1680  // Sine part P_k^S has no contribution
1681 
1682  } // End of loop over pressure shape functions
1683  } // End of Jacobian calculation
1684 
1685  } // End of if not boundary condition statement
1686 
1687  } // End of loop over shape functions
1688 
1689 
1690  // ====================
1691  // CONTINUITY EQUATIONS
1692  // ====================
1693 
1694  // Loop over the pressure shape functions
1695  for (unsigned l = 0; l < n_pres; l++)
1696  {
1697  // Cache test function
1698  const double testp_ = testp[l];
1699 
1700  // --------------------------------------
1701  // FIRST CONTINUITY EQUATION: COSINE PART
1702  // --------------------------------------
1703 
1704  // Get local equation number of first pressure value at this node
1705  local_eqn = p_local_eqn(l, 0);
1706 
1707  // If it's not a boundary condition
1708  if (local_eqn >= 0)
1709  {
1710  // Gradient terms
1711  residuals[local_eqn] +=
1712  (interpolated_UC + r * interpolated_dUCdr + k * interpolated_VS +
1713  r * interpolated_dWCdz) *
1714  testp_ * W;
1715 
1716  // Calculate the Jacobian
1717  // ----------------------
1718 
1719  if (flag)
1720  {
1721  // Loop over the velocity shape functions
1722  for (unsigned l2 = 0; l2 < n_node; l2++)
1723  {
1724  // Cache velocity shape functions and their derivatives
1725  const double psif_ = psif[l2];
1726  const double dpsifdr = dpsifdx(l2, 0);
1727  const double dpsifdz = dpsifdx(l2, 1);
1728 
1729  // Radial velocity component (cosine part) U_k^C
1730  local_unknown = nodal_local_eqn(l2, u_nodal_index[0]);
1731  if (local_unknown >= 0)
1732  {
1733  jacobian(local_eqn, local_unknown) +=
1734  (psif_ + r * dpsifdr) * testp_ * W;
1735  }
1736 
1737  // Radial velocity component (sine part) U_k^S
1738  // has no contribution
1739 
1740  // Axial velocity component (cosine part) W_k^C
1741  local_unknown = nodal_local_eqn(l2, u_nodal_index[2]);
1742  if (local_unknown >= 0)
1743  {
1744  jacobian(local_eqn, local_unknown) += r * dpsifdz * testp_ * W;
1745  }
1746 
1747  // Axial velocity component (sine part) W_k^S
1748  // has no contribution
1749 
1750  // Azimuthal velocity component (cosine part) V_k^C
1751  // has no contribution
1752 
1753  // Azimuthal velocity component (sine part) V_k^S
1754  local_unknown = nodal_local_eqn(l2, u_nodal_index[5]);
1755  if (local_unknown >= 0)
1756  {
1757  jacobian(local_eqn, local_unknown) += k * psif_ * testp_ * W;
1758  }
1759  } // End of loop over velocity shape functions
1760 
1761  // Real and imaginary pressure components, P_k^C and P_k^S,
1762  // have no contribution
1763 
1764  } // End of Jacobian calculation
1765 
1766  } // End of if not boundary condition statement
1767 
1768  // -------------------------------------
1769  // SECOND CONTINUITY EQUATION: SINE PART
1770  // -------------------------------------
1771 
1772  // Get local equation number of second pressure value at this node
1773  local_eqn = p_local_eqn(l, 1);
1774 
1775  // If it's not a boundary condition
1776  if (local_eqn >= 0)
1777  {
1778  // Gradient terms
1779  residuals[local_eqn] +=
1780  (interpolated_US + r * interpolated_dUSdr - k * interpolated_VC +
1781  r * interpolated_dWSdz) *
1782  testp_ * W;
1783 
1784  // Calculate the Jacobian
1785  // ----------------------
1786 
1787  if (flag)
1788  {
1789  // Loop over the velocity shape functions
1790  for (unsigned l2 = 0; l2 < n_node; l2++)
1791  {
1792  // Cache velocity shape functions and their derivatives
1793  const double psif_ = psif[l2];
1794  const double dpsifdr = dpsifdx(l2, 0);
1795  const double dpsifdz = dpsifdx(l2, 1);
1796 
1797  // Radial velocity component (cosine part) U_k^C
1798  // has no contribution
1799 
1800  // Radial velocity component (sine part) U_k^S
1801  local_unknown = nodal_local_eqn(l2, u_nodal_index[1]);
1802  if (local_unknown >= 0)
1803  {
1804  jacobian(local_eqn, local_unknown) +=
1805  (psif_ + r * dpsifdr) * testp_ * W;
1806  }
1807 
1808  // Axial velocity component (cosine part) W_k^C
1809  // has no contribution
1810 
1811  // Axial velocity component (sine part) W_k^S
1812  local_unknown = nodal_local_eqn(l2, u_nodal_index[3]);
1813  if (local_unknown >= 0)
1814  {
1815  jacobian(local_eqn, local_unknown) += r * dpsifdz * testp_ * W;
1816  }
1817 
1818  // Azimuthal velocity component (cosine part) V_k^C
1819  local_unknown = nodal_local_eqn(l2, u_nodal_index[4]);
1820  if (local_unknown >= 0)
1821  {
1822  jacobian(local_eqn, local_unknown) -= k * psif_ * testp_ * W;
1823  }
1824 
1825  // Azimuthal velocity component (sine part) V_k^S
1826  // has no contribution
1827 
1828  } // End of loop over velocity shape functions
1829 
1830  // Real and imaginary pressure components, P_k^C and P_k^S,
1831  // have no contribution
1832 
1833  } // End of Jacobian calculation
1834 
1835  } // End of if not boundary condition statement
1836 
1837  } // End of loop over pressure shape functions
1838 
1839  } // End of loop over the integration points
1840 
1841  } // End of fill_in_generic_residual_contribution_linearised_axi_nst
1842 
1843 
1844  /// ////////////////////////////////////////////////////////////////////////
1845  /// ////////////////////////////////////////////////////////////////////////
1846  /// ////////////////////////////////////////////////////////////////////////
1847 
1848 
1849  /// Linearised axisymmetric Crouzeix-Raviart elements
1850  /// -------------------------------------------------
1851 
1852 
1853  //=======================================================================
1854  /// Set the data for the number of variables at each node
1855  //=======================================================================
1856  const unsigned
1858  6, 6, 6, 6, 6, 6, 6, 6, 6};
1859 
1860 
1861  //========================================================================
1862  /// Number of values (pinned or dofs) required at node n
1863  //========================================================================
1865  const unsigned& n) const
1866  {
1867  return Initial_Nvalue[n];
1868  }
1869 
1870 
1871  /// ////////////////////////////////////////////////////////////////////////
1872  /// ////////////////////////////////////////////////////////////////////////
1873  /// ////////////////////////////////////////////////////////////////////////
1874 
1875 
1876  /// Linearised axisymmetric Taylor-Hood elements
1877  /// --------------------------------------------
1878 
1879 
1880  //=======================================================================
1881  /// Set the data for the number of variables at each node
1882  //=======================================================================
1884  8, 6, 8, 6, 6, 6, 8, 6, 8};
1885 
1886 
1887  //=======================================================================
1888  /// Set the data for the pressure conversion array
1889  //=======================================================================
1891  0, 2, 6, 8};
1892 
1893 
1894 } // 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 long nrow() const
Return the number of rows of the matrix.
Definition: matrices.h:485
unsigned long ncol() const
Return the number of columns of the matrix.
Definition: matrices.h:491
Node *& node_pt(const unsigned &n)
Return a pointer to the local node n.
Definition: elements.h:2179
double nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n. Produces suitably interpolated values for hanging nodes...
Definition: elements.h:2597
virtual std::string tecplot_zone_string(const unsigned &nplot) const
Return string for tecplot zone header (when plotting nplot points in each "coordinate direction")
Definition: elements.h:3165
virtual 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:3992
int nodal_local_eqn(const unsigned &n, const unsigned &i) const
Return the local equation number corresponding to the i-th value at the n-th local node.
Definition: elements.h:1436
unsigned nnode() const
Return the number of nodes.
Definition: elements.h:2214
Integral *const & integral_pt() const
Return the pointer to the integration scheme (const version)
Definition: elements.h:1967
virtual void get_s_plot(const unsigned &i, const unsigned &nplot, Vector< double > &s, const bool &shifted_to_interior=false) const
Get cector of local coordinates of plot point i (when plotting nplot points in each "coordinate direc...
Definition: elements.h:3152
virtual unsigned nplot_points(const unsigned &nplot) const
Return total number of plot points (when plotting nplot points in each "coordinate direction")
Definition: elements.h:3190
double 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:2321
double dshape_eulerian(const Vector< double > &s, Shape &psi, DShape &dpsidx) const
Compute the geometric shape functions and also first derivatives w.r.t. global coordinates at local c...
Definition: elements.cc:3328
double raw_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:2260
double raw_nodal_value(const unsigned &n, const unsigned &i) const
Return the i-th value stored at local node n but do NOT take hanging nodes into account.
Definition: elements.h:2580
double raw_nodal_position(const unsigned &n, const unsigned &i) const
Return the i-th coordinate at local node n. Do not use the hanging node representation....
Definition: elements.cc:1714
virtual void write_tecplot_zone_footer(std::ostream &outfile, const unsigned &nplot) const
Add tecplot zone "footer" to output stream (when plotting nplot points in each "coordinate direction"...
Definition: elements.h:3178
virtual 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.
const double & viscosity_ratio() const
Viscosity ratio for element: element's viscosity relative to the viscosity used in the definition of ...
const double & re_st() const
Product of Reynolds and Strouhal number (=Womersley number)
static Vector< double > Gamma
Vector to decide whether the stress-divergence form is used or not.
double du_dt_linearised_axi_nst(const unsigned &n, const unsigned &i) const
Return the i-th component of du/dt at local node n. Uses suitably interpolated value for hanging node...
void output(std::ostream &outfile)
Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, P^C, P^S in tecplot format. Default number of pl...
static double Default_Physical_Ratio_Value
Static default value for the physical ratios (all initialised to one)
virtual void get_base_flow_u(const double &time, const unsigned &ipt, const Vector< double > &x, Vector< double > &result) const
Calculate the velocity components of the base flow solution at a given time and Eulerian position.
static int Pressure_not_stored_at_node
Static "magic" number that indicates that the pressure is not stored at a node.
virtual int p_local_eqn(const unsigned &n, const unsigned &i)=0
Access function for the local equation number information for the i-th component of the pressure....
bool ALE_is_disabled
Boolean flag to indicate if ALE formulation is disabled when the time-derivatives are computed....
virtual void pshape_linearised_axi_nst(const Vector< double > &s, Shape &psi) const =0
Compute the pressure shape functions at local coordinate s.
double interpolated_u_linearised_axi_nst(const Vector< double > &s, const unsigned &i) const
Return the i-th component of the FE interpolated velocity u[i] at local coordinate s.
void output_veloc(std::ostream &outfile, const unsigned &nplot, const unsigned &t)
Output function: r, z, U^C, U^S, V^C, V^S, W^C, W^S, in tecplot format. nplot points in each coordina...
virtual void get_base_flow_dudx(const double &time, const unsigned &ipt, const Vector< double > &x, DenseMatrix< double > &result) const
Calculate the derivatives of the velocity components of the base flow solution w.r....
virtual double dshape_and_dtest_eulerian_at_knot_linearised_axi_nst(const unsigned &ipt, Shape &psi, DShape &dpsidx, Shape &test, DShape &dtestdx) const =0
Compute the shape functions and their derivatives w.r.t. global coordinates at the ipt-th integration...
virtual unsigned npres_linearised_axi_nst() const =0
Return the number of pressure degrees of freedom associated with a single pressure component in the e...
double interpolated_p_linearised_axi_nst(const Vector< double > &s, const unsigned &i) const
Return the i-th component of the FE interpolated pressure p[i] at local coordinate s.
const double & density_ratio() const
Density ratio for element: element's density relative to the viscosity used in the definition of the ...
virtual unsigned u_index_linearised_axi_nst(const unsigned &i) const
Return the index at which the i-th unknown velocity component is stored. The default value,...
const int & azimuthal_mode_number() const
Azimuthal mode number k in e^ik(theta) decomposition.
virtual void fill_in_generic_residual_contribution_linearised_axi_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...
virtual double p_linearised_axi_nst(const unsigned &n_p, const unsigned &i) const =0
Return the i-th pressure value at local pressure "node" n_p. Uses suitably interpolated value for han...
void strain_rate(const Vector< double > &s, DenseMatrix< double > &strain_rate) const
Strain-rate tensor: where (in that order)
virtual unsigned required_nvalue(const unsigned &n) const
Return number of values (pinned or dofs) required at local node n.
static const unsigned Initial_Nvalue[]
Static array of ints to hold required number of variables at nodes.
static const unsigned Initial_Nvalue[]
Static array of ints to hold number of variables at node.
static const unsigned Pconv[]
Static array of ints to hold conversion from pressure node numbers to actual node numbers.
An OomphLibError object which should be thrown when an run-time error is encountered....
A Class for shape functions. In simple cases, the shape functions have only one index that can be tho...
Definition: shape.h:76
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
double & time()
Return current value of continous time.
Definition: timesteppers.h:332
//////////////////////////////////////////////////////////////////// ////////////////////////////////...